ALEPH_1996_S3196992.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/RivetYODA.hh" 00004 #include "Rivet/Tools/Logging.hh" 00005 #include "Rivet/Projections/FinalState.hh" 00006 #include "Rivet/Projections/ChargedFinalState.hh" 00007 #include "Rivet/Projections/IdentifiedFinalState.hh" 00008 #include "Rivet/Projections/FastJets.hh" 00009 #include "Rivet/Projections/Thrust.hh" 00010 00011 namespace Rivet { 00012 00013 00014 /// @brief ALEPH measurement of quark-to-photon fragmentation function 00015 class ALEPH_1996_S3196992 : public Analysis { 00016 public: 00017 00018 /// Constructor 00019 ALEPH_1996_S3196992() : Analysis("ALEPH_1996_S3196992") 00020 { 00021 } 00022 00023 00024 /// @name Analysis methods 00025 //@{ 00026 00027 void init() { 00028 // Set up projections 00029 FinalState fs; 00030 addProjection(FastJets(fs, FastJets::DURHAM, 0.7), "DurhamJets"); 00031 IdentifiedFinalState ifs(-MAXRAPIDITY, +MAXRAPIDITY, 0.0); 00032 ifs.acceptId(PHOTON); 00033 addProjection(ifs, "Photons"); 00034 addProjection(Thrust(fs), "Thrust"); 00035 addProjection(ChargedFinalState(), "CFS"); 00036 00037 // Book histograms 00038 _h_z_2jet_001 = bookHisto1D(1, 1, 1); 00039 _h_z_2jet_006 = bookHisto1D(2, 1, 1); 00040 _h_z_2jet_01 = bookHisto1D(3, 1, 1); 00041 _h_z_2jet_033 = bookHisto1D(4, 1, 1); 00042 _h_z_3jet_001 = bookHisto1D(5, 1, 1); 00043 _h_z_3jet_006 = bookHisto1D(6, 1, 1); 00044 _h_z_3jet_01 = bookHisto1D(7, 1, 1); 00045 _h_z_4jet_001 = bookHisto1D(8, 1, 1); 00046 } 00047 00048 00049 /// Perform the per-event analysis 00050 void analyze(const Event& event) { 00051 const double weight = event.weight(); 00052 00053 if (applyProjection<FinalState>(event, "CFS").particles().size()<2) { 00054 vetoEvent; 00055 } 00056 00057 const ParticleVector allphotons = applyProjection<IdentifiedFinalState>(event, "Photons").particles(); 00058 ParticleVector photons; 00059 foreach (const Particle& photon, allphotons) { 00060 if (fabs(cos(photon.momentum().theta()))<0.95 && photon.momentum().E()>5.0*GeV) { 00061 photons.push_back(photon); 00062 } 00063 } 00064 if (photons.size()<1) { 00065 vetoEvent; 00066 } 00067 00068 const Thrust& thrust = applyProjection<Thrust>(event, "Thrust"); 00069 if (fabs(cos(thrust.thrustAxis().theta()))>0.95) { 00070 vetoEvent; 00071 } 00072 00073 const FastJets& durjet = applyProjection<FastJets>(event, "DurhamJets"); 00074 00075 foreach (const Particle& photon, photons) { 00076 00077 PseudoJets jets_001 = durjet.clusterSeq()->exclusive_jets_ycut(0.01); 00078 foreach (const fastjet::PseudoJet& jet, jets_001) { 00079 if (particleInJet(photon, jet, durjet.clusterSeq())) { 00080 double zgamma = photon.momentum().E()/jet.E(); 00081 if (jets_001.size() == 2) _h_z_2jet_001->fill(zgamma, weight); 00082 else if (jets_001.size() == 3) _h_z_3jet_001->fill(zgamma, weight); 00083 else if (jets_001.size() > 3) _h_z_4jet_001->fill(zgamma, weight); 00084 break; 00085 } 00086 } 00087 00088 PseudoJets jets_006 = durjet.clusterSeq()->exclusive_jets_ycut(0.06); 00089 foreach (const fastjet::PseudoJet& jet, jets_006) { 00090 if (particleInJet(photon, jet, durjet.clusterSeq())) { 00091 double zgamma = photon.momentum().E()/jet.E(); 00092 if (jets_006.size() == 2) _h_z_2jet_006->fill(zgamma, weight); 00093 else if (jets_006.size() == 3) _h_z_3jet_006->fill(zgamma, weight); 00094 break; 00095 } 00096 } 00097 00098 PseudoJets jets_01 = durjet.clusterSeq()->exclusive_jets_ycut(0.1); 00099 foreach (const fastjet::PseudoJet& jet, jets_01) { 00100 if (particleInJet(photon, jet, durjet.clusterSeq())) { 00101 double zgamma = photon.momentum().E()/jet.E(); 00102 if (jets_01.size() == 2) _h_z_2jet_01->fill(zgamma, weight); 00103 else if (jets_01.size() == 3) _h_z_3jet_01->fill(zgamma, weight); 00104 break; 00105 } 00106 } 00107 00108 PseudoJets jets_033 = durjet.clusterSeq()->exclusive_jets_ycut(0.33); 00109 foreach (const fastjet::PseudoJet& jet, jets_033) { 00110 if (particleInJet(photon, jet, durjet.clusterSeq())) { 00111 double zgamma = photon.momentum().E()/jet.E(); 00112 if (jets_033.size() == 2) _h_z_2jet_033->fill(zgamma, weight); 00113 break; 00114 } 00115 } 00116 00117 } 00118 } 00119 00120 00121 bool particleInJet(const Particle& p, const fastjet::PseudoJet& jet, 00122 const fastjet::ClusterSequence* cseq ) { 00123 foreach (const fastjet::PseudoJet& jetpart, cseq->constituents(jet)) { 00124 if (fuzzyEquals(jetpart.E(), p.momentum().E()) && 00125 fuzzyEquals(jetpart.px(), p.momentum().x()) && 00126 fuzzyEquals(jetpart.py(), p.momentum().y()) && 00127 fuzzyEquals(jetpart.pz(), p.momentum().z())) { 00128 return true; 00129 } 00130 } 00131 return false; 00132 } 00133 00134 00135 /// Normalise histograms etc., after the run 00136 void finalize() { 00137 scale(_h_z_2jet_001, 1000.0/sumOfWeights()); 00138 scale(_h_z_2jet_006, 1000.0/sumOfWeights()); 00139 scale(_h_z_2jet_01, 1000.0/sumOfWeights()); 00140 scale(_h_z_2jet_033, 1000.0/sumOfWeights()); 00141 scale(_h_z_3jet_001, 1000.0/sumOfWeights()); 00142 scale(_h_z_3jet_006, 1000.0/sumOfWeights()); 00143 scale(_h_z_3jet_01, 1000.0/sumOfWeights()); 00144 scale(_h_z_4jet_001, 1000.0/sumOfWeights()); 00145 } 00146 00147 //@} 00148 00149 00150 private: 00151 00152 /// @name Histograms 00153 //@{ 00154 00155 Histo1DPtr _h_z_2jet_001, _h_z_2jet_006, _h_z_2jet_01, _h_z_2jet_033; 00156 Histo1DPtr _h_z_3jet_001, _h_z_3jet_006, _h_z_3jet_01; 00157 Histo1DPtr _h_z_4jet_001; 00158 //@} 00159 00160 }; 00161 00162 00163 00164 // The hook for the plugin system 00165 DECLARE_RIVET_PLUGIN(ALEPH_1996_S3196992); 00166 00167 } Generated on Fri Dec 21 2012 15:03:38 for The Rivet MC analysis system by ![]() |