00001
00002 #include "Rivet/Analyses/MC_JetAnalysis.hh"
00003 #include "Rivet/Projections/FastJets.hh"
00004 #include "Rivet/Tools/Logging.hh"
00005 #include "Rivet/RivetAIDA.hh"
00006
00007 namespace Rivet {
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), m_njet(njet), m_jetpro_name(jetpro_name), m_jetptcut(jetptcut),
00015 _h_log10_d(njet, NULL), _h_log10_R(njet+1, NULL), _h_pT_jet(njet, NULL),
00016 _h_eta_jet(njet, NULL), _h_eta_jet_plus(njet), _h_eta_jet_minus(njet),
00017 _h_rap_jet(njet, NULL), _h_rap_jet_plus(njet), _h_rap_jet_minus(njet),
00018 _h_mass_jet(njet, NULL)
00019 {
00020 setNeedsCrossSection(true);
00021 }
00022
00023
00024
00025
00026 void MC_JetAnalysis::init() {
00027
00028 for (size_t i=0; i < m_njet; ++i) {
00029 stringstream dname;
00030 dname << "log10_d_" << i << i+1;
00031
00032 _h_log10_d[i] = bookHistogram1D(dname.str(), 50, 0.2, log10(0.5*sqrtS()));
00033
00034 stringstream Rname;
00035 Rname << "log10_R_" << i;
00036 _h_log10_R[i] = bookDataPointSet(Rname.str(), 50, 0.2, log10(0.5*sqrtS()));
00037
00038 stringstream pTname;
00039 pTname << "jet_pT_" << i+1;
00040 double pTmax = 1.0/(double(i)+2.0) * sqrtS()/GeV/2.0;
00041 int nbins_pT = 100/(i+1);
00042 _h_pT_jet[i] = bookHistogram1D(pTname.str(), logBinEdges(nbins_pT, 10.0, pTmax));
00043
00044 stringstream massname;
00045 massname << "jet_mass_" << i+1;
00046 double mmax = 100.0;
00047 int nbins_m = 100/(i+1);
00048 _h_mass_jet[i] = bookHistogram1D(massname.str(), logBinEdges(nbins_m, 1.0, mmax));
00049
00050 stringstream etaname;
00051 etaname << "jet_eta_" << i+1;
00052 _h_eta_jet[i] = bookHistogram1D(etaname.str(), i > 1 ? 25 : 50, -5.0, 5.0);
00053 _h_eta_jet_plus[i].reset(new LWH::Histogram1D(i > 1 ? 15 : 25, 0, 5));
00054 _h_eta_jet_minus[i].reset(new LWH::Histogram1D(i > 1 ? 15 : 25, 0, 5));
00055
00056 stringstream rapname;
00057 rapname << "jet_y_" << i+1;
00058 _h_rap_jet[i] = bookHistogram1D(rapname.str(), i>1 ? 25 : 50, -5.0, 5.0);
00059 _h_rap_jet_plus[i].reset(new LWH::Histogram1D(i > 1 ? 15 : 25, 0, 5));
00060 _h_rap_jet_minus[i].reset(new LWH::Histogram1D(i > 1 ? 15 : 25, 0, 5));
00061
00062 for (size_t j = i+1; j < m_njet; ++j) {
00063 std::pair<size_t, size_t> ij = std::make_pair(i, j);
00064
00065 stringstream detaname;
00066 detaname << "jets_deta_" << i+1 << j+1;
00067 _h_deta_jets.insert(make_pair(ij, bookHistogram1D(detaname.str(), 25, -5.0, 5.0)));
00068
00069 stringstream dRname;
00070 dRname << "jets_dR_" << i+1 << j+1;
00071 _h_dR_jets.insert(make_pair(ij, bookHistogram1D(dRname.str(), 25, 0.0, 5.0)));
00072 }
00073 }
00074 stringstream Rname;
00075 Rname << "log10_R_" << m_njet;
00076 _h_log10_R[m_njet] = bookDataPointSet(Rname.str(), 50, 0.2, log10(0.5*sqrtS()));
00077
00078 _h_jet_multi_exclusive = bookHistogram1D("jet_multi_exclusive", m_njet+3, -0.5, m_njet+3-0.5);
00079 _h_jet_multi_inclusive = bookHistogram1D("jet_multi_inclusive", m_njet+3, -0.5, m_njet+3-0.5);
00080 _h_jet_multi_ratio = bookDataPointSet("jet_multi_ratio", m_njet+2, 0.5, m_njet+3-0.5);
00081 _h_jet_HT = bookHistogram1D("jet_HT", logBinEdges(50, m_jetptcut, sqrtS()/GeV/2.0));
00082 }
00083
00084
00085
00086
00087 void MC_JetAnalysis::analyze(const Event & e) {
00088 const double weight = e.weight();
00089
00090 const FastJets& jetpro = applyProjection<FastJets>(e, m_jetpro_name);
00091
00092
00093 const fastjet::ClusterSequence* seq = jetpro.clusterSeq();
00094 if (seq != NULL) {
00095 double previous_dij = 10.0;
00096 for (size_t i = 0; i < m_njet; ++i) {
00097
00098 double d_ij = log10(sqrt(seq->exclusive_dmerge_max(i)));
00099
00100
00101 _h_log10_d[i]->fill(d_ij, weight);
00102
00103
00104 for (int ibin = 0; ibin < _h_log10_R[i]->size(); ++ibin) {
00105 IDataPoint* dp = _h_log10_R[i]->point(ibin);
00106 double dcut = dp->coordinate(0)->value();
00107 if (d_ij < dcut && previous_dij > dcut) {
00108 dp->coordinate(1)->setValue(dp->coordinate(1)->value() + weight);
00109 }
00110 }
00111 previous_dij = d_ij;
00112 }
00113
00114 for (int ibin = 0; ibin<_h_log10_R[m_njet]->size(); ++ibin) {
00115 IDataPoint* dp = _h_log10_R[m_njet]->point(ibin);
00116 double dcut = dp->coordinate(0)->value();
00117 if (previous_dij > dcut) {
00118 dp->coordinate(1)->setValue(dp->coordinate(1)->value() + weight);
00119 }
00120 }
00121 }
00122
00123 const Jets& jets = jetpro.jetsByPt(m_jetptcut);
00124
00125
00126 for (size_t i = 0; i < m_njet; ++i) {
00127 if (jets.size() < i+1) continue;
00128 _h_pT_jet[i]->fill(jets[i].momentum().pT()/GeV, weight);
00129
00130 double m2_i = jets[i].momentum().mass2();
00131 if (m2_i < 0) {
00132 if (m2_i < -1e-4) {
00133 getLog() << Log::WARNING << "Jet mass2 is negative: " << m2_i << " GeV^2. "
00134 << "Truncating to 0.0, assuming numerical precision is to blame." << endl;
00135 }
00136 m2_i = 0.0;
00137 }
00138
00139
00140 _h_mass_jet[i]->fill(sqrt(m2_i)/GeV, weight);
00141
00142
00143 const double eta_i = jets[i].momentum().eta();
00144 _h_eta_jet[i]->fill(eta_i, weight);
00145 if (eta_i > 0.0) {
00146 _h_eta_jet_plus[i]->fill(eta_i, weight);
00147 } else {
00148 _h_eta_jet_minus[i]->fill(fabs(eta_i), weight);
00149 }
00150
00151
00152 const double rap_i = jets[i].momentum().rapidity();
00153 _h_rap_jet[i]->fill(rap_i, weight);
00154 if (rap_i > 0.0) {
00155 _h_rap_jet_plus[i]->fill(rap_i, weight);
00156 } else {
00157 _h_rap_jet_minus[i]->fill(fabs(rap_i), weight);
00158 }
00159
00160
00161 for (size_t j = i+1; j < m_njet; ++j) {
00162 if (jets.size() < j+1) continue;
00163 std::pair<size_t, size_t> ij = std::make_pair(i, j);
00164 double deta = jets[i].momentum().eta()-jets[j].momentum().eta();
00165 double dR = deltaR(jets[i].momentum(), jets[j].momentum());
00166 _h_deta_jets[ij]->fill(deta, weight);
00167 _h_dR_jets[ij]->fill(dR, weight);
00168 }
00169 }
00170 _h_jet_multi_exclusive->fill(jets.size(), weight);
00171
00172 for (size_t i = 0; i < m_njet+2; ++i) {
00173 if (jets.size() >= i) {
00174 _h_jet_multi_inclusive->fill(i, weight);
00175 }
00176 }
00177
00178 double HT=0.0;
00179 foreach (const Jet& jet, jets) {
00180 HT += jet.momentum().pT();
00181 }
00182 _h_jet_HT->fill(HT, weight);
00183 }
00184
00185
00186
00187 void MC_JetAnalysis::finalize() {
00188 for (size_t i = 0; i < m_njet; ++i) {
00189 scale(_h_log10_d[i], crossSection()/sumOfWeights());
00190 for (int ibin = 0; ibin<_h_log10_R[i]->size(); ++ibin) {
00191 IDataPoint* dp = _h_log10_R[i]->point(ibin);
00192 dp->coordinate(1)->setValue(dp->coordinate(1)->value()*crossSection()/sumOfWeights());
00193 }
00194
00195 scale(_h_pT_jet[i], crossSection()/sumOfWeights());
00196 scale(_h_mass_jet[i], crossSection()/sumOfWeights());
00197 scale(_h_eta_jet[i], crossSection()/sumOfWeights());
00198 scale(_h_rap_jet[i], crossSection()/sumOfWeights());
00199
00200
00201 stringstream etaname;
00202 etaname << "jet_eta_pmratio_" << i+1;
00203 histogramFactory().divide(histoPath(etaname.str()), *_h_eta_jet_plus[i], *_h_eta_jet_minus[i]);
00204 stringstream rapname;
00205 rapname << "jet_y_pmratio_" << i+1;
00206 histogramFactory().divide(histoPath(rapname.str()), *_h_rap_jet_plus[i], *_h_rap_jet_minus[i]);
00207 }
00208
00209 for (int ibin = 0; ibin < _h_log10_R[m_njet]->size(); ++ibin) {
00210 IDataPoint* dp =_h_log10_R[m_njet]->point(ibin);
00211 dp->coordinate(1)->setValue(dp->coordinate(1)->value()*crossSection()/sumOfWeights());
00212 }
00213
00214
00215 map<pair<size_t, size_t>, AIDA::IHistogram1D*>::iterator it;
00216 for (it = _h_deta_jets.begin(); it != _h_deta_jets.end(); ++it) {
00217 scale(it->second, crossSection()/sumOfWeights());
00218 }
00219 for (it = _h_dR_jets.begin(); it != _h_dR_jets.end(); ++it) {
00220 scale(it->second, crossSection()/sumOfWeights());
00221 }
00222
00223
00224 int Nbins = _h_jet_multi_inclusive->axis().bins();
00225 std::vector<double> ratio(Nbins-1, 0.0);
00226 std::vector<double> err(Nbins-1, 0.0);
00227 for (int i = 0; i < Nbins-1; ++i) {
00228 if (_h_jet_multi_inclusive->binHeight(i) > 0.0 && _h_jet_multi_inclusive->binHeight(i+1) > 0.0) {
00229 ratio[i] = _h_jet_multi_inclusive->binHeight(i+1)/_h_jet_multi_inclusive->binHeight(i);
00230 double relerr_i = _h_jet_multi_inclusive->binError(i)/_h_jet_multi_inclusive->binHeight(i);
00231 double relerr_j = _h_jet_multi_inclusive->binError(i+1)/_h_jet_multi_inclusive->binHeight(i+1);
00232 err[i] = ratio[i] * (relerr_i + relerr_j);
00233 }
00234 }
00235 _h_jet_multi_ratio->setCoordinate(1, ratio, err);
00236
00237 scale(_h_jet_multi_exclusive, crossSection()/sumOfWeights());
00238 scale(_h_jet_multi_inclusive, crossSection()/sumOfWeights());
00239 scale(_h_jet_HT, crossSection()/sumOfWeights());
00240 }
00241
00242
00243 }