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