ATLAS_2010_S8894728.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/RivetYODA.hh" 00004 #include "Rivet/Projections/ChargedFinalState.hh" 00005 #include "Rivet/Tools/Logging.hh" 00006 00007 namespace Rivet { 00008 00009 00010 class ATLAS_2010_S8894728 : public Analysis { 00011 public: 00012 00013 ATLAS_2010_S8894728() : Analysis("ATLAS_2010_S8894728") { } 00014 00015 00016 void init() { 00017 const ChargedFinalState cfs100(-2.5, 2.5, 100*MeV); 00018 addProjection(cfs100, "CFS100"); 00019 const ChargedFinalState cfs500(-2.5, 2.5, 500*MeV); 00020 addProjection(cfs500, "CFS500"); 00021 const ChargedFinalState cfslead(-2.5, 2.5, 1.0*GeV); 00022 addProjection(cfslead, "CFSlead"); 00023 00024 // Get an index for the beam energy 00025 int isqrts = -1; 00026 if (fuzzyEquals(sqrtS(), 900*GeV)) isqrts = 0; 00027 else if (fuzzyEquals(sqrtS(), 7*TeV)) isqrts = 1; 00028 assert(isqrts >= 0); 00029 00030 // Nch profiles, 500 MeV track pT cut 00031 _hist_nch_transverse_500[0] = bookProfile1D(1+isqrts, 1, 1); 00032 _hist_nch_toward_500 = bookProfile1D(1+isqrts, 1, 2); 00033 _hist_nch_away_500 = bookProfile1D(1+isqrts, 1, 3); 00034 00035 // pTsum profiles, 500 MeV track pT cut 00036 _hist_ptsum_transverse_500[0] = bookProfile1D(3+isqrts, 1, 1); 00037 _hist_ptsum_toward_500 = bookProfile1D(3+isqrts, 1, 2); 00038 _hist_ptsum_away_500 = bookProfile1D(3+isqrts, 1, 3); 00039 00040 // Standard deviation profiles 00041 // First the higher moments of main profiles to calculate variance and error on variance... 00042 for (size_t i = 1; i < 4; ++i) { 00043 _hist_nch_transverse_500[i].reset(new Profile1D(referenceData(1+isqrts, 1, 1))); 00044 _hist_ptsum_transverse_500[i].reset(new Profile1D(referenceData(3+isqrts, 1, 1))); 00045 } 00046 // Then the data point sets into which the results will be inserted 00047 _dps_sdnch_transverse_500 = bookScatter2D(5+isqrts, 1, 1); 00048 _dps_sdptsum_transverse_500 = bookScatter2D(7+isqrts, 1, 1); 00049 00050 // <pT> profiles, 500 MeV track pT cut 00051 _hist_ptavg_transverse_500 = bookProfile1D(9+isqrts, 1, 1); 00052 _hist_ptavg_toward_500 = bookProfile1D(9+isqrts, 1, 2); 00053 _hist_ptavg_away_500 = bookProfile1D(9+isqrts, 1, 3); 00054 00055 // <pT> vs. Nch profiles, 500 MeV track pT cut 00056 _hist_dn_dpt_transverse_500 = bookProfile1D(11+isqrts, 1, 1); 00057 _hist_dn_dpt_toward_500 = bookProfile1D(11+isqrts, 1, 2); 00058 _hist_dn_dpt_away_500 = bookProfile1D(11+isqrts, 1, 3); 00059 00060 // Nch vs. Delta(phi) profiles, 500 MeV track pT cut 00061 _hist_N_vs_dPhi_1_500 = bookProfile1D(13+isqrts, 1, 1); 00062 _hist_N_vs_dPhi_2_500 = bookProfile1D(13+isqrts, 1, 2); 00063 _hist_N_vs_dPhi_3_500 = bookProfile1D(13+isqrts, 1, 3); 00064 _hist_N_vs_dPhi_5_500 = bookProfile1D(13+isqrts, 1, 4); 00065 // pT vs. Delta(phi) profiles, 500 MeV track pT cut 00066 _hist_pT_vs_dPhi_1_500 = bookProfile1D(15+isqrts, 1, 1); 00067 _hist_pT_vs_dPhi_2_500 = bookProfile1D(15+isqrts, 1, 2); 00068 _hist_pT_vs_dPhi_3_500 = bookProfile1D(15+isqrts, 1, 3); 00069 _hist_pT_vs_dPhi_5_500 = bookProfile1D(15+isqrts, 1, 4); 00070 00071 // Nch and pTsum profiles, 100 MeV track pT cut 00072 _hist_nch_transverse_100 = bookProfile1D(17+isqrts, 1, 1); 00073 _hist_nch_toward_100 = bookProfile1D(17+isqrts, 1, 2); 00074 _hist_nch_away_100 = bookProfile1D(17+isqrts, 1, 3); 00075 _hist_ptsum_transverse_100 = bookProfile1D(19+isqrts, 1, 1); 00076 _hist_ptsum_toward_100 = bookProfile1D(19+isqrts, 1, 2); 00077 _hist_ptsum_away_100 = bookProfile1D(19+isqrts, 1, 3); 00078 00079 // Profiles vs. eta (7 TeV only) 00080 if (isqrts == 1) { 00081 _hist_nch_vs_eta_transverse_100 = bookProfile1D(21, 1, 1); 00082 _hist_ptsum_vs_eta_transverse_100 = bookProfile1D(22, 1, 1); 00083 } 00084 00085 } 00086 00087 00088 // Little helper function to identify Delta(phi) regions 00089 inline int region_index(double dphi) { 00090 assert(inRange(dphi, 0.0, PI, CLOSED, CLOSED)); 00091 if (dphi < PI/3.0) return 0; 00092 if (dphi < 2*PI/3.0) return 1; 00093 return 2; 00094 } 00095 00096 00097 void analyze(const Event& event) { 00098 const double weight = event.weight(); 00099 00100 // Require at least one track in the event with pT >= 1 GeV 00101 const ChargedFinalState& cfslead = applyProjection<ChargedFinalState>(event, "CFSlead"); 00102 if (cfslead.size() < 1) { 00103 vetoEvent; 00104 } 00105 00106 // These are the charged particles (tracks) with pT > 500 MeV 00107 const ChargedFinalState& charged500 = applyProjection<ChargedFinalState>(event, "CFS500"); 00108 00109 // Identify leading track and its phi and pT (this is the same for both the 100 MeV and 500 MeV track cuts) 00110 ParticleVector particles500 = charged500.particlesByPt(); 00111 Particle p_lead = particles500[0]; 00112 const double philead = p_lead.momentum().phi(); 00113 const double etalead = p_lead.momentum().eta(); 00114 const double pTlead = p_lead.momentum().perp(); 00115 MSG_DEBUG("Leading track: pT = " << pTlead << ", eta = " << etalead << ", phi = " << philead); 00116 00117 // Iterate over all > 500 MeV particles and count particles and scalar pTsum in the three regions 00118 vector<double> num500(3, 0), ptSum500(3, 0.0); 00119 // Temporary histos that bin Nch and pT in dPhi. 00120 // NB. Only one of each needed since binnings are the same for the energies and pT cuts 00121 Histo1D hist_num_dphi_500(referenceData(13,1,1)); 00122 Histo1D hist_pt_dphi_500(referenceData(15,1,1)); 00123 foreach (const Particle& p, particles500) { 00124 const double pT = p.momentum().pT(); 00125 const double dPhi = deltaPhi(philead, p.momentum().phi()); 00126 const int ir = region_index(dPhi); 00127 num500[ir] += 1; 00128 ptSum500[ir] += pT; 00129 00130 // Fill temp histos to bin Nch and pT in dPhi 00131 if (p.genParticle() != p_lead.genParticle()) { // We don't want to fill all those zeros from the leading track... 00132 hist_num_dphi_500.fill(dPhi, 1); 00133 hist_pt_dphi_500.fill(dPhi, pT); 00134 } 00135 } 00136 00137 00138 // Iterate over charged particles again for profiles against Nch 00139 // This is necessary since the Nch are region-specific and so are only known after the first loop 00140 foreach (const Particle& p, particles500) { 00141 const double pT = p.momentum().pT(); 00142 const double dPhi = deltaPhi(philead, p.momentum().phi()); 00143 const int ir = region_index(dPhi); 00144 switch (ir) { 00145 case 0: 00146 _hist_dn_dpt_toward_500->fill(num500[0], pT, weight); 00147 break; 00148 case 1: 00149 _hist_dn_dpt_transverse_500->fill(num500[1], pT, weight); 00150 break; 00151 case 2: 00152 _hist_dn_dpt_away_500->fill(num500[2], pT, weight); 00153 break; 00154 default: 00155 assert(false && "How did we get here?"); 00156 } 00157 } 00158 00159 00160 // Now fill underlying event histograms 00161 // The densities are calculated by dividing the UE properties by dEta*dPhi 00162 // -- each region has a dPhi of 2*PI/3 and dEta is two times 2.5 00163 const double dEtadPhi = (2*2.5 * 2*PI/3.0); 00164 // Transverse profiles need 4 orders of moments for stddev with errors 00165 for (int i = 0; i < 4; ++i) { 00166 _hist_nch_transverse_500[i]->fill(pTlead/GeV, intpow(num500[1]/dEtadPhi, i+1), weight); 00167 _hist_ptsum_transverse_500[i]->fill(pTlead/GeV, intpow(ptSum500[1]/GeV/dEtadPhi, i+1), weight); 00168 } 00169 // Toward and away profiles only need the first moment (the mean) 00170 _hist_nch_toward_500->fill(pTlead/GeV, num500[0]/dEtadPhi, weight); 00171 _hist_nch_away_500->fill(pTlead/GeV, num500[2]/dEtadPhi, weight); 00172 _hist_ptsum_toward_500->fill(pTlead/GeV, ptSum500[0]/GeV/dEtadPhi, weight); 00173 _hist_ptsum_away_500->fill(pTlead/GeV, ptSum500[2]/GeV/dEtadPhi, weight); 00174 // <pT> profiles 00175 //MSG_INFO("Trans pT1, pTsum, Nch, <pT>" << pTlead/GeV << ", " << ptSum500[1]/GeV << ", " << num500[1] << ", " << ptSum500[1]/GeV/num500[1]); 00176 if (num500[1] > 0) _hist_ptavg_transverse_500->fill(pTlead/GeV, ptSum500[1]/GeV/num500[1], weight); 00177 if (num500[0] > 0) _hist_ptavg_toward_500->fill(pTlead/GeV, ptSum500[0]/GeV/num500[0], weight); 00178 if (num500[2] > 0) _hist_ptavg_away_500->fill(pTlead/GeV, ptSum500[2]/GeV/num500[2], weight); 00179 00180 00181 // Update the "proper" dphi profile histograms 00182 // Note that we fill dN/dEtadPhi: dEta = 2*2.5, dPhi = 2*PI/nBins 00183 // The values tabulated in the note are for an (undefined) signed Delta(phi) rather than 00184 // |Delta(phi)| and so differ by a factor of 2: we have to actually norm for angular range = 2pi 00185 const size_t nbins = referenceData(13,1,1).numPoints(); 00186 std::vector<double> ptcut; 00187 if (fuzzyEquals(sqrtS(), 900*GeV)) { 00188 ptcut += 1.0; ptcut += 1.5; ptcut += 2.0; ptcut += 2.5; 00189 } 00190 else if (fuzzyEquals(sqrtS(), 7*TeV)) { 00191 ptcut += 1.0; ptcut += 2.0; ptcut += 3.0; ptcut += 5.0; 00192 } 00193 assert(ptcut.size() == 4); 00194 for (size_t i = 0; i < nbins; ++i) { 00195 // First Nch 00196 double mean = hist_num_dphi_500.bin(i).midpoint(); 00197 double value = 0.; 00198 if (hist_num_dphi_500.bin(i).numEntries() > 0) { 00199 mean = hist_num_dphi_500.bin(i).xMean(); 00200 value = hist_num_dphi_500.bin(i).area()/hist_num_dphi_500.bin(i).width()/10.0; 00201 } 00202 if (pTlead/GeV >= ptcut[0]) _hist_N_vs_dPhi_1_500->fill(mean, value, weight); 00203 if (pTlead/GeV >= ptcut[1]) _hist_N_vs_dPhi_2_500->fill(mean, value, weight); 00204 if (pTlead/GeV >= ptcut[2]) _hist_N_vs_dPhi_3_500->fill(mean, value, weight); 00205 if (pTlead/GeV >= ptcut[3]) _hist_N_vs_dPhi_5_500->fill(mean, value, weight); 00206 00207 // Then pT 00208 mean = hist_pt_dphi_500.bin(i).midpoint(); 00209 value = 0.; 00210 if (hist_pt_dphi_500.bin(i).numEntries() > 0) { 00211 mean = hist_pt_dphi_500.bin(i).xMean(); 00212 value = hist_pt_dphi_500.bin(i).area()/hist_pt_dphi_500.bin(i).width()/10.0; 00213 } 00214 if (pTlead/GeV >= ptcut[0]) _hist_pT_vs_dPhi_1_500->fill(mean, value, weight); 00215 if (pTlead/GeV >= ptcut[1]) _hist_pT_vs_dPhi_2_500->fill(mean, value, weight); 00216 if (pTlead/GeV >= ptcut[2]) _hist_pT_vs_dPhi_3_500->fill(mean, value, weight); 00217 if (pTlead/GeV >= ptcut[3]) _hist_pT_vs_dPhi_5_500->fill(mean, value, weight); 00218 } 00219 00220 00221 ////////////////////// 00222 00223 00224 // These are the charged particles (tracks) with pT > 100 MeV 00225 const ChargedFinalState& charged100 = applyProjection<ChargedFinalState>(event, "CFS100"); 00226 00227 // Iterate over all > 100 MeV particles and count particles and scalar pTsum in the three regions 00228 vector<double> num100(3, 0), ptSum100(3, 0.0); 00229 foreach (const Particle& p, charged100.particles()) { 00230 const double pT = p.momentum().pT(); 00231 const double dPhi = deltaPhi(philead, p.momentum().phi()); 00232 const int ir = region_index(dPhi); 00233 num100[ir] += 1; 00234 ptSum100[ir] += pT; 00235 } 00236 00237 // Now fill the two sets of 100 MeV underlying event histograms 00238 _hist_nch_transverse_100->fill(pTlead/GeV, num100[1]/dEtadPhi, weight); 00239 _hist_nch_toward_100->fill(pTlead/GeV, num100[0]/dEtadPhi, weight); 00240 _hist_nch_away_100->fill(pTlead/GeV, num100[2]/dEtadPhi, weight); 00241 _hist_ptsum_transverse_100->fill(pTlead/GeV, ptSum100[1]/GeV/dEtadPhi, weight); 00242 _hist_ptsum_toward_100->fill(pTlead/GeV, ptSum100[0]/GeV/dEtadPhi, weight); 00243 _hist_ptsum_away_100->fill(pTlead/GeV, ptSum100[2]/GeV/dEtadPhi, weight); 00244 00245 // And finally the Nch and pT vs eta_lead profiles (again from > 100 MeV tracks, and only at 7 TeV) 00246 if (fuzzyEquals(sqrtS(), 7*TeV) && pTlead > 5*GeV) { 00247 // MSG_INFO(sqrtS() << " " << pTlead << " " << ptSum100[1]/dEtadPhi << " " << num100[1]/dEtadPhi); 00248 _hist_nch_vs_eta_transverse_100->fill(etalead, num100[1]/dEtadPhi, weight); 00249 _hist_ptsum_vs_eta_transverse_100->fill(etalead, ptSum100[1]/GeV/dEtadPhi, weight); 00250 } 00251 00252 } 00253 00254 00255 void finalize() { 00256 // Convert the various moments of the 500 MeV trans pT and Nch distributions to std devs with correct error 00257 _moments_to_stddev(_hist_nch_transverse_500, _dps_sdnch_transverse_500); 00258 _moments_to_stddev(_hist_ptsum_transverse_500, _dps_sdptsum_transverse_500); 00259 } 00260 00261 00262 private: 00263 00264 00265 inline void _moments_to_stddev(Profile1DPtr moment_profiles[], Scatter2DPtr target_dps) { 00266 for (size_t b = 0; b < moment_profiles[0]->numBins(); ++b) { // loop over points 00267 /// @todo Assuming unit weights here! Should use N_effective = sumW**2/sumW2? 00268 const double numentries = moment_profiles[0]->bin(b).numEntries(); 00269 const double x = moment_profiles[0]->bin(b).midpoint(); 00270 const double ex = moment_profiles[0]->bin(b).width()/2.; 00271 double var = 0.; 00272 double sd = 0.; 00273 if (numentries > 0) { 00274 var = moment_profiles[1]->bin(b).mean() - intpow(moment_profiles[0]->bin(b).mean(), 2); 00275 sd = fuzzyLessEquals(var,0.) ? 0 : sqrt(var); //< Numerical safety check 00276 } 00277 if (sd == 0 || numentries < 3) { 00278 MSG_WARNING("Need at least 3 bin entries and a non-zero central value to calculate " 00279 << "an error on standard deviation profiles (bin " << b << ")"); 00280 target_dps->addPoint(x, sd, ex, 0); 00281 continue; 00282 } 00283 // c2(y) = m4(x) - 4 m3(x) m1(x) - m2(x)^2 + 8 m2(x) m1(x)^2 - 4 m1(x)^4 00284 const double var_on_var = moment_profiles[3]->bin(b).mean() 00285 - 4 * moment_profiles[2]->bin(b).mean() * moment_profiles[0]->bin(b).mean() 00286 - intpow(moment_profiles[1]->bin(b).mean(), 2) 00287 + 8 * moment_profiles[1]->bin(b).mean() * intpow(moment_profiles[0]->bin(b).mean(), 2) 00288 - 4 * intpow(moment_profiles[0]->bin(b).mean(), 4); 00289 const double stderr_on_var = sqrt(var_on_var/(numentries-2.0)); 00290 const double stderr_on_sd = stderr_on_var / (2.0*sd); 00291 target_dps->addPoint(x, sd, ex, stderr_on_sd); 00292 } 00293 } 00294 00295 00296 private: 00297 00298 Profile1DPtr _hist_nch_transverse_500[4]; 00299 Profile1DPtr _hist_nch_toward_500; 00300 Profile1DPtr _hist_nch_away_500; 00301 00302 Profile1DPtr _hist_ptsum_transverse_500[4]; 00303 Profile1DPtr _hist_ptsum_toward_500; 00304 Profile1DPtr _hist_ptsum_away_500; 00305 00306 Scatter2DPtr _dps_sdnch_transverse_500; 00307 Scatter2DPtr _dps_sdptsum_transverse_500; 00308 00309 Profile1DPtr _hist_ptavg_transverse_500; 00310 Profile1DPtr _hist_ptavg_toward_500; 00311 Profile1DPtr _hist_ptavg_away_500; 00312 00313 Profile1DPtr _hist_dn_dpt_transverse_500; 00314 Profile1DPtr _hist_dn_dpt_toward_500; 00315 Profile1DPtr _hist_dn_dpt_away_500; 00316 00317 Profile1DPtr _hist_N_vs_dPhi_1_500; 00318 Profile1DPtr _hist_N_vs_dPhi_2_500; 00319 Profile1DPtr _hist_N_vs_dPhi_3_500; 00320 Profile1DPtr _hist_N_vs_dPhi_5_500; 00321 00322 Profile1DPtr _hist_pT_vs_dPhi_1_500; 00323 Profile1DPtr _hist_pT_vs_dPhi_2_500; 00324 Profile1DPtr _hist_pT_vs_dPhi_3_500; 00325 Profile1DPtr _hist_pT_vs_dPhi_5_500; 00326 00327 Profile1DPtr _hist_nch_transverse_100; 00328 Profile1DPtr _hist_nch_toward_100; 00329 Profile1DPtr _hist_nch_away_100; 00330 Profile1DPtr _hist_ptsum_transverse_100; 00331 Profile1DPtr _hist_ptsum_toward_100; 00332 Profile1DPtr _hist_ptsum_away_100; 00333 00334 Profile1DPtr _hist_nch_vs_eta_transverse_100; 00335 Profile1DPtr _hist_ptsum_vs_eta_transverse_100; 00336 00337 }; 00338 00339 00340 00341 // The hook for the plugin system 00342 DECLARE_RIVET_PLUGIN(ATLAS_2010_S8894728); 00343 00344 } Generated on Fri Dec 21 2012 15:03:38 for The Rivet MC analysis system by ![]() |