CMS_2013_I1261026.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/FinalState.hh" 00004 #include "Rivet/Projections/ChargedFinalState.hh" 00005 #include "Rivet/Projections/FastJets.hh" 00006 #include "Rivet/Projections/Beam.hh" 00007 #include "Rivet/Projections/VetoedFinalState.hh" 00008 00009 namespace Rivet { 00010 00011 00012 /// Jet and underlying event properties as a function of particle multiplicity 00013 class CMS_2013_I1261026 : public Analysis { 00014 public: 00015 00016 CMS_2013_I1261026() 00017 : Analysis("CMS_2013_I1261026"), _jetStructNorm(5,0.), _multBinCent(5,0.), 00018 _jetCounter5GeV(5,0.), _jetCounter30GeV(5,0.), _passedEv(5,0.) 00019 { } 00020 00021 00022 void init() { 00023 FastJets jetpro(ChargedFinalState(-2.4, 2.4, 0.25*GeV), FastJets::ANTIKT, 0.5); 00024 addProjection(jetpro, "Jets"); 00025 00026 const ChargedFinalState cfs(-2.4, 2.4, 0.25*GeV); 00027 addProjection(cfs, "CFS250"); 00028 00029 // For min bias trigger 00030 const ChargedFinalState cfsBSCplus(3.23, 4.65, 500*MeV); 00031 addProjection(cfsBSCplus, "cfsBSCplus"); 00032 00033 const ChargedFinalState cfsBSCminus(-4.65, -3.23, 500*MeV); 00034 addProjection(cfsBSCminus, "cfsBSCminus"); 00035 00036 // Histograms: 00037 _h_AllTrkMeanPt = bookProfile1D(1, 1, 1); 00038 _h_SoftTrkMeanPt = bookProfile1D(2, 1, 1); 00039 _h_IntrajetTrkMeanPt = bookProfile1D(3, 1, 1); 00040 _h_IntrajetLeaderTrkMeanPt = bookProfile1D(4, 1, 1); 00041 _h_MeanJetPt = bookProfile1D(5, 1, 1); 00042 _h_JetRate5GeV = bookProfile1D(6, 1, 1); 00043 _h_JetRate30GeV = bookProfile1D(7, 1, 1); 00044 00045 for (int ihist = 0; ihist < 5; ++ihist) { 00046 _h_JetSpectrum[ihist] = bookHisto1D(ihist+8, 1, 1); 00047 _h_JetStruct[ihist] = bookHisto1D(ihist+13, 1, 1); 00048 00049 // Temp histograms for distribution parameters and SEM calculation 00050 _th_AllTrkSpectrum[ihist] = Histo1D(200, 0.0, 20.0); 00051 _th_SoftTrkSpectrum[ihist] = Histo1D(100, 0.0, 15.0); 00052 _th_JetTrkSpectrum[ihist] = Histo1D(100, 0.0, 20.0); 00053 _th_JetLTrkSpectrum[ihist] = Histo1D(100, 0.0, 20.0); 00054 } 00055 00056 _multBinCent[0] = 20; 00057 _multBinCent[1] = 40; 00058 _multBinCent[2] = 65; 00059 _multBinCent[3] = 95; 00060 _multBinCent[4] = 125; 00061 } 00062 00063 00064 /// Perform the per-event analysis 00065 void analyze(const Event& event) { 00066 const double weight = event.weight(); 00067 00068 // MinBias trigger 00069 const ChargedFinalState& cfsBSCplus = applyProjection<ChargedFinalState>(event, "cfsBSCplus"); 00070 if (cfsBSCplus.empty()) vetoEvent; 00071 const ChargedFinalState& cfsBSCminus = applyProjection<ChargedFinalState>(event, "cfsBSCminus"); 00072 if (cfsBSCminus.empty()) vetoEvent; 00073 00074 const ChargedFinalState& cfsp = applyProjection<ChargedFinalState>(event, "CFS250"); 00075 if (cfsp.empty()) vetoEvent; 00076 00077 const FastJets& jetpro = applyProjection<FastJets>(event, "Jets"); 00078 const Jets& jets = jetpro.jetsByPt(5.0*GeV); 00079 00080 const int mult = cfsp.size(); 00081 00082 int multbin[6] = { 10, 30, 50, 80, 110, 140 }; 00083 for (int ibin = 0; ibin < 5; ++ibin) { 00084 if (mult > multbin[ibin] && mult <= multbin[ibin + 1]) { 00085 _passedEv[ibin] += weight; 00086 eventDecomp(event, ibin, weight); 00087 00088 for (size_t ijets = 0; ijets < jets.size(); ijets++) { 00089 if (jets[ijets].abseta() < 1.9) { 00090 _h_JetSpectrum[ibin]->fill(jets[ijets].pT()/GeV, weight); 00091 if (jets[ijets].pT() > 5*GeV) _jetCounter5GeV[ibin] += weight; 00092 if (jets[ijets].pT() > 30*GeV) _jetCounter30GeV[ibin] += weight; 00093 } 00094 } 00095 } 00096 } 00097 00098 } 00099 00100 00101 /// Normalise histograms etc., after the run 00102 void finalize() { 00103 for (size_t i = 0; i < 5; ++i) { 00104 // All trk mean pT vs Nch 00105 _h_AllTrkMeanPt->fill(_multBinCent[i], _th_AllTrkSpectrum[i].xMean(), getMeanError(_th_AllTrkSpectrum[i])); 00106 00107 // Soft trk mean pT vs Nch 00108 _h_SoftTrkMeanPt->fill(_multBinCent[i], _th_SoftTrkSpectrum[i].xMean(), getMeanError(_th_SoftTrkSpectrum[i])); 00109 00110 // Intrajet trk mean pT vs Nch 00111 _h_IntrajetTrkMeanPt->fill(_multBinCent[i], _th_JetTrkSpectrum[i].xMean(), getMeanError(_th_JetTrkSpectrum[i])); 00112 00113 // Intrajet leader trk mean pT vs Nch 00114 _h_IntrajetLeaderTrkMeanPt->fill(_multBinCent[i], _th_JetLTrkSpectrum[i].xMean(), getMeanError(_th_JetLTrkSpectrum[i])); 00115 00116 // Jet mean pT vs Nch 00117 const double sem = (_h_JetSpectrum[i]->xStdDev())/(sqrt(_h_JetSpectrum[i]->sumW())) / _h_JetSpectrum[i]->xMean(); 00118 _h_MeanJetPt->fill(_multBinCent[i], _h_JetSpectrum[i]->xMean(), sem); 00119 00120 // Jet rates 00121 double avJetRate5 = _jetCounter5GeV[i] / _passedEv[i]; 00122 double avJetRate30 = _jetCounter30GeV[i] / _passedEv[i]; 00123 00124 const double sem5 = (_jetCounter5GeV[i] != 0) ? 1 / sqrt(_jetCounter5GeV[i]) : 0; 00125 _h_JetRate5GeV->fill(_multBinCent[i], avJetRate5, sem5); 00126 00127 const double sem30 = (_jetCounter30GeV[i] != 0) ? 1 / sqrt(_jetCounter30GeV[i]) : 0; 00128 _h_JetRate30GeV->fill(_multBinCent[i], avJetRate30, sem30); 00129 00130 scale(_h_JetSpectrum[i], 4.0 / _jetCounter5GeV[i]); 00131 scale(_h_JetStruct[i], 0.08 / _jetStructNorm[i]); 00132 } 00133 00134 } 00135 00136 00137 double getMeanError(const Histo1D& hist) { 00138 double sem = hist.xStdErr(); // Standard error of the mean 00139 return sem / hist.xMean(); // relative SEM 00140 } 00141 00142 00143 void eventDecomp(const Event& event, size_t ibin, double weight) { 00144 00145 struct TrkInJet { double pt; double eta; double phi; double R; }; 00146 TrkInJet jetConstituents[100][100]; //1-st index - the number of the jet, 2-nd index - track in the jet 00147 TrkInJet jetsEv[100]; 00148 size_t j[100]; 00149 size_t jCount = 0; 00150 00151 for (size_t i = 0; i < 100; i++) { 00152 j[i] = 0; 00153 jetsEv[i].pt = 0; 00154 jetsEv[i].eta = 0; 00155 jetsEv[i].phi = 0; 00156 for (size_t k = 0; k < 100; k++) { 00157 jetConstituents[i][k].pt = 0; 00158 jetConstituents[i][k].phi = 0; 00159 jetConstituents[i][k].eta = 0; 00160 jetConstituents[i][k].R = 0; 00161 } 00162 } 00163 00164 const FastJets& jetpro = applyProjection<FastJets>(event, "Jets"); 00165 const Jets& jets = jetpro.jetsByPt(5.0*GeV); 00166 00167 // Start event decomp 00168 00169 for (size_t ijets = 0; ijets < jets.size(); ijets++) { 00170 jetsEv[ijets].pt = jets[ijets].pT(); 00171 jetsEv[ijets].eta = jets[ijets].eta(); 00172 jetsEv[ijets].phi = jets[ijets].phi(); 00173 jCount++; 00174 } 00175 00176 const ChargedFinalState& cfsp = applyProjection<ChargedFinalState>(event, "CFS250"); 00177 foreach (const Particle& p, cfsp.particles()) { 00178 _th_AllTrkSpectrum[ibin].fill(p.pT()/GeV, weight); 00179 int flag = 0; 00180 for (size_t i = 0; i < jCount; i++) { 00181 const double delta_phi = deltaPhi(jetsEv[i].phi, p.phi()); 00182 const double delta_eta = jetsEv[i].eta - p.eta(); 00183 const double R = sqrt(delta_phi * delta_phi + delta_eta * delta_eta); 00184 if (R <= 0.5) { 00185 flag++; 00186 jetConstituents[i][j[i]].pt = p.pT(); 00187 jetConstituents[i][j[i]].R = R; 00188 j[i]++; 00189 } 00190 } 00191 if (flag == 0) _th_SoftTrkSpectrum[ibin].fill(p.pT()/GeV, weight); 00192 } 00193 00194 for (size_t i = 0; i < jCount; i++) { 00195 double ptInjetLeader = 0; 00196 if (!inRange(jetsEv[i].eta, -1.9, 1.9)) continue; // only fully reconstructed jets for internal jet studies 00197 for (size_t k = 0; k < j[i]; k++) { 00198 _th_JetTrkSpectrum[ibin].fill(jetConstituents[i][k].pt , weight); 00199 _h_JetStruct[ibin]->fill(jetConstituents[i][k].R, jetConstituents[i][k].pt/jetsEv[i].pt); 00200 _jetStructNorm[ibin] += jetConstituents[i][k].pt / jetsEv[i].pt; 00201 if (ptInjetLeader < jetConstituents[i][k].pt) ptInjetLeader = jetConstituents[i][k].pt; 00202 } 00203 if (ptInjetLeader != 0) _th_JetLTrkSpectrum[ibin].fill(ptInjetLeader, weight); 00204 } 00205 00206 } 00207 00208 00209 private: 00210 00211 // Counters etc. 00212 vector<double> _jetStructNorm; 00213 vector<double> _multBinCent; 00214 /// @todo Need to handle weights 00215 vector<double> _jetCounter5GeV, _jetCounter30GeV, _passedEv; 00216 00217 Profile1DPtr _h_AllTrkMeanPt, _h_SoftTrkMeanPt; 00218 Profile1DPtr _h_IntrajetTrkMeanPt, _h_IntrajetLeaderTrkMeanPt; 00219 Profile1DPtr _h_MeanJetPt; 00220 Profile1DPtr _h_JetRate5GeV, _h_JetRate30GeV; 00221 00222 Histo1DPtr _h_JetSpectrum[5]; 00223 Histo1DPtr _h_JetStruct[5]; 00224 00225 // Temp histograms 00226 Histo1D _th_AllTrkSpectrum[5]; 00227 Histo1D _th_SoftTrkSpectrum[5]; 00228 Histo1D _th_JetTrkSpectrum[5]; 00229 Histo1D _th_JetLTrkSpectrum[5]; 00230 00231 }; 00232 00233 00234 // The hook for the plugin system 00235 DECLARE_RIVET_PLUGIN(CMS_2013_I1261026); 00236 00237 } Generated on Wed Oct 7 2015 12:09:12 for The Rivet MC analysis system by ![]() |