MC_JetAnalysis.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analyses/MC_JetAnalysis.hh" 00003 #include "Rivet/Projections/FastJets.hh" 00004 00005 namespace Rivet { 00006 00007 00008 00009 00010 MC_JetAnalysis::MC_JetAnalysis(const string& name, 00011 size_t njet, 00012 const string& jetpro_name, 00013 double jetptcut) 00014 : Analysis(name), _njet(njet), _jetpro_name(jetpro_name), _jetptcut(jetptcut), 00015 _h_pT_jet(njet), 00016 _h_eta_jet(njet), _h_eta_jet_plus(njet), _h_eta_jet_minus(njet), 00017 _h_rap_jet(njet), _h_rap_jet_plus(njet), _h_rap_jet_minus(njet), 00018 _h_mass_jet(njet) 00019 { 00020 setNeedsCrossSection(true); // legitimate use, since a base class has no .info file! 00021 } 00022 00023 00024 00025 // Book histograms 00026 void MC_JetAnalysis::init() { 00027 const double sqrts = sqrtS() ? sqrtS() : 14000.*GeV; 00028 00029 for (size_t i = 0; i < _njet; ++i) { 00030 const string pTname = "jet_pT_" + to_str(i+1); 00031 const double pTmax = 1.0/(double(i)+2.0) * sqrts/GeV/2.0; 00032 const int nbins_pT = 100/(i+1); 00033 _h_pT_jet[i] = bookHisto1D(pTname, logspace(nbins_pT, 10.0, pTmax)); 00034 00035 const string massname = "jet_mass_" + to_str(i+1); 00036 const double mmax = 100.0; 00037 const int nbins_m = 100/(i+1); 00038 _h_mass_jet[i] = bookHisto1D(massname, logspace(nbins_m, 1.0, mmax)); 00039 00040 const string etaname = "jet_eta_" + to_str(i+1); 00041 _h_eta_jet[i] = bookHisto1D(etaname, i > 1 ? 25 : 50, -5.0, 5.0); 00042 _h_eta_jet_plus[i].reset(new Histo1D(i > 1 ? 15 : 25, 0, 5)); 00043 _h_eta_jet_minus[i].reset(new Histo1D(i > 1 ? 15 : 25, 0, 5)); 00044 00045 const string rapname = "jet_y_" + to_str(i+1); 00046 _h_rap_jet[i] = bookHisto1D(rapname, i>1 ? 25 : 50, -5.0, 5.0); 00047 _h_rap_jet_plus[i].reset(new Histo1D(i > 1 ? 15 : 25, 0, 5)); 00048 _h_rap_jet_minus[i].reset(new Histo1D(i > 1 ? 15 : 25, 0, 5)); 00049 00050 for (size_t j = i+1; j < min(size_t(3), _njet); ++j) { 00051 const std::pair<size_t, size_t> ij = std::make_pair(i, j); 00052 00053 string detaname = "jets_deta_" + to_str(i+1) + to_str(j+1); 00054 _h_deta_jets.insert(make_pair(ij, bookHisto1D(detaname, 25, -5.0, 5.0))); 00055 00056 string dphiname = "jets_dphi_" + to_str(i+1) + to_str(j+1); 00057 _h_dphi_jets.insert(make_pair(ij, bookHisto1D(dphiname, 25, 0.0, M_PI))); 00058 00059 string dRname = "jets_dR_" + to_str(i+1) + to_str(j+1); 00060 _h_dR_jets.insert(make_pair(ij, bookHisto1D(dRname, 25, 0.0, 5.0))); 00061 } 00062 } 00063 00064 _h_jet_multi_exclusive = bookHisto1D("jet_multi_exclusive", _njet+3, -0.5, _njet+3-0.5); 00065 _h_jet_multi_inclusive = bookHisto1D("jet_multi_inclusive", _njet+3, -0.5, _njet+3-0.5); 00066 _h_jet_multi_ratio = bookScatter2D("jet_multi_ratio"); 00067 _h_jet_HT = bookHisto1D("jet_HT", logspace(50, _jetptcut, sqrts/GeV/2.0)); 00068 } 00069 00070 00071 // Do the analysis 00072 void MC_JetAnalysis::analyze(const Event & e) { 00073 const double weight = e.weight(); 00074 00075 const Jets& jets = applyProjection<FastJets>(e, _jetpro_name).jetsByPt(_jetptcut); 00076 00077 for (size_t i = 0; i < _njet; ++i) { 00078 if (jets.size() < i+1) continue; 00079 _h_pT_jet[i]->fill(jets[i].pT()/GeV, weight); 00080 // Check for numerical precision issues with jet masses 00081 double m2_i = jets[i].mass2(); 00082 if (m2_i < 0) { 00083 if (m2_i < -1e-4) { 00084 MSG_WARNING("Jet mass2 is negative: " << m2_i << " GeV^2. " 00085 << "Truncating to 0.0, assuming numerical precision is to blame."); 00086 } 00087 m2_i = 0.0; 00088 } 00089 00090 // Jet mass 00091 _h_mass_jet[i]->fill(sqrt(m2_i)/GeV, weight); 00092 00093 // Jet eta 00094 const double eta_i = jets[i].eta(); 00095 _h_eta_jet[i]->fill(eta_i, weight); 00096 (eta_i > 0.0 ? _h_eta_jet_plus : _h_eta_jet_minus)[i]->fill(fabs(eta_i), weight); 00097 00098 // Jet rapidity 00099 const double rap_i = jets[i].rapidity(); 00100 _h_rap_jet[i]->fill(rap_i, weight); 00101 (rap_i > 0.0 ? _h_rap_jet_plus : _h_rap_jet_minus)[i]->fill(fabs(rap_i), weight); 00102 00103 // Inter-jet properties 00104 for (size_t j = i+1; j < min(size_t(3),_njet); ++j) { 00105 if (jets.size() < j+1) continue; 00106 std::pair<size_t, size_t> ij = std::make_pair(i, j); 00107 double deta = jets[i].eta()-jets[j].eta(); 00108 double dphi = deltaPhi(jets[i].momentum(),jets[j].momentum()); 00109 double dR = deltaR(jets[i].momentum(), jets[j].momentum()); 00110 _h_deta_jets[ij]->fill(deta, weight); 00111 _h_dphi_jets[ij]->fill(dphi, weight); 00112 _h_dR_jets[ij]->fill(dR, weight); 00113 } 00114 } 00115 00116 // Multiplicities 00117 _h_jet_multi_exclusive->fill(jets.size(), weight); 00118 for (size_t i = 0; i < _njet+2; ++i) { 00119 if (jets.size() >= i) { 00120 _h_jet_multi_inclusive->fill(i, weight); 00121 } 00122 } 00123 00124 // HT 00125 double HT = 0.0; 00126 foreach (const Jet& jet, jets) { 00127 HT += jet.pT(); 00128 } 00129 _h_jet_HT->fill(HT, weight); 00130 } 00131 00132 00133 // Finalize 00134 void MC_JetAnalysis::finalize() { 00135 for (size_t i = 0; i < _njet; ++i) { 00136 scale(_h_pT_jet[i], crossSection()/sumOfWeights()); 00137 scale(_h_mass_jet[i], crossSection()/sumOfWeights()); 00138 scale(_h_eta_jet[i], crossSection()/sumOfWeights()); 00139 scale(_h_rap_jet[i], crossSection()/sumOfWeights()); 00140 00141 // Create eta/rapidity ratio plots 00142 divide(*_h_eta_jet_plus[i], *_h_eta_jet_minus[i], bookScatter2D("jet_eta_pmratio_" + to_str(i+1))); 00143 divide(*_h_rap_jet_plus[i], *_h_rap_jet_minus[i], bookScatter2D("jet_y_pmratio_" + to_str(i+1))); 00144 } 00145 00146 // Scale the d{eta,phi,R} histograms 00147 typedef map<pair<size_t, size_t>, Histo1DPtr> HistMap; 00148 foreach (HistMap::value_type& it, _h_deta_jets) scale(it.second, crossSection()/sumOfWeights()); 00149 foreach (HistMap::value_type& it, _h_dphi_jets) scale(it.second, crossSection()/sumOfWeights()); 00150 foreach (HistMap::value_type& it, _h_dR_jets) scale(it.second, crossSection()/sumOfWeights()); 00151 00152 // Fill inclusive jet multi ratio 00153 int Nbins = _h_jet_multi_inclusive->numBins(); 00154 for (int i = 0; i < Nbins-1; ++i) { 00155 _h_jet_multi_ratio->addPoint(i+1, 0, 0.5, 0); 00156 if (_h_jet_multi_inclusive->bin(i).sumW() > 0.0) { 00157 const double ratio = _h_jet_multi_inclusive->bin(i+1).sumW()/_h_jet_multi_inclusive->bin(i).sumW(); 00158 const double relerr_i = _h_jet_multi_inclusive->bin(i).relErr(); 00159 const double relerr_j = _h_jet_multi_inclusive->bin(i+1).relErr(); 00160 const double err = ratio * (relerr_i + relerr_j); 00161 _h_jet_multi_ratio->point(i).setY(ratio, err); 00162 } 00163 } 00164 00165 scale(_h_jet_multi_exclusive, crossSection()/sumOfWeights()); 00166 scale(_h_jet_multi_inclusive, crossSection()/sumOfWeights()); 00167 scale(_h_jet_HT, crossSection()/sumOfWeights()); 00168 } 00169 00170 00171 } Generated on Thu Mar 10 2016 08:29:51 for The Rivet MC analysis system by ![]() |