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 #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   // Book histograms
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   }
00082 
00083 
00084 
00085   // Do the analysis
00086   void MC_JetAnalysis::analyze(const Event & e) {
00087     const double weight = e.weight();
00088 
00089     const FastJets& jetpro = applyProjection<FastJets>(e, m_jetpro_name);
00090 
00091     // Jet resolutions and integrated jet rates
00092     const fastjet::ClusterSequence* seq = jetpro.clusterSeq();
00093     if (seq != NULL) {
00094       double previous_dij = 10.0;
00095       for (size_t i = 0; i < m_njet; ++i) {
00096         // Jet resolution i -> j
00097         double d_ij = log10(sqrt(seq->exclusive_dmerge_max(i)));
00098 
00099         // Fill differential jet resolution
00100         _h_log10_d[i]->fill(d_ij, weight);
00101 
00102         // Fill integrated jet resolution
00103         for (int ibin = 0; ibin < _h_log10_R[i]->size(); ++ibin) {
00104           IDataPoint* dp = _h_log10_R[i]->point(ibin);
00105           double dcut = dp->coordinate(0)->value();
00106           if (d_ij < dcut && previous_dij > dcut) {
00107             dp->coordinate(1)->setValue(dp->coordinate(1)->value() + weight);
00108           }
00109         }
00110         previous_dij = d_ij;
00111       }
00112       // One remaining integrated jet resolution
00113       for (int ibin = 0; ibin<_h_log10_R[m_njet]->size(); ++ibin) {
00114         IDataPoint* dp = _h_log10_R[m_njet]->point(ibin);
00115         double dcut = dp->coordinate(0)->value();
00116         if (previous_dij > dcut) {
00117           dp->coordinate(1)->setValue(dp->coordinate(1)->value() + weight);
00118         }
00119       }
00120     }
00121 
00122     const Jets& jets = jetpro.jetsByPt(m_jetptcut);
00123 
00124     // The remaining direct jet observables
00125     for (size_t i = 0; i < m_njet; ++i) {
00126       if (jets.size() < i+1) continue;
00127       _h_pT_jet[i]->fill(jets[i].momentum().pT()/GeV, weight);
00128       // Check for numerical precision issues with jet masses
00129       double m2_i = jets[i].momentum().mass2();
00130       if (m2_i < 0) {
00131         if (m2_i < -1e-8) {
00132           getLog() << Log::WARNING << "Jet mass2 is negative: " << m2_i << " GeV^2. "
00133                    << "Truncating to 0.0, assuming numerical precision is to blame." << endl;
00134         }
00135         m2_i = 0.0;
00136       }
00137 
00138       // Jet mass
00139       _h_mass_jet[i]->fill(sqrt(m2_i)/GeV, weight);
00140 
00141       // Jet eta
00142       const double eta_i = jets[i].momentum().eta();
00143       _h_eta_jet[i]->fill(eta_i, weight);
00144       if (eta_i > 0.0) {
00145         _h_eta_jet_plus[i]->fill(eta_i, weight);
00146       } else {
00147         _h_eta_jet_minus[i]->fill(fabs(eta_i), weight);
00148       }
00149 
00150       // Jet rapidity
00151       const double rap_i = jets[i].momentum().rapidity();
00152       _h_rap_jet[i]->fill(rap_i, weight);
00153       if (rap_i > 0.0) {
00154         _h_rap_jet_plus[i]->fill(rap_i, weight);
00155       } else {
00156         _h_rap_jet_minus[i]->fill(fabs(rap_i), weight);
00157       }
00158 
00159       // Inter-jet properties
00160       for (size_t j = i+1; j < m_njet; ++j) {
00161         if (jets.size() < j+1) continue;
00162         std::pair<size_t, size_t> ij = std::make_pair(i, j);
00163         double deta = jets[i].momentum().eta()-jets[j].momentum().eta();
00164         double dR = deltaR(jets[i].momentum(), jets[j].momentum());
00165         _h_deta_jets[ij]->fill(deta, weight);
00166         _h_dR_jets[ij]->fill(dR, weight);
00167       }
00168     }
00169     _h_jet_multi_exclusive->fill(jets.size(), weight);
00170 
00171     for (size_t i = 0; i < m_njet+2; ++i) {
00172       if (jets.size() >= i) {
00173         _h_jet_multi_inclusive->fill(i, weight);
00174       }
00175     }
00176   }
00177 
00178 
00179   // Finalize
00180   void MC_JetAnalysis::finalize() {
00181     for (size_t i = 0; i < m_njet; ++i) {
00182       scale(_h_log10_d[i], crossSection()/sumOfWeights());
00183       for (int ibin = 0; ibin<_h_log10_R[i]->size(); ++ibin) {
00184         IDataPoint* dp = _h_log10_R[i]->point(ibin);
00185         dp->coordinate(1)->setValue(dp->coordinate(1)->value()*crossSection()/sumOfWeights());
00186       }
00187 
00188       scale(_h_pT_jet[i], crossSection()/sumOfWeights());
00189       scale(_h_mass_jet[i], crossSection()/sumOfWeights());
00190       scale(_h_eta_jet[i], crossSection()/sumOfWeights());
00191       scale(_h_rap_jet[i], crossSection()/sumOfWeights());
00192 
00193       // Create eta/rapidity ratio plots
00194       stringstream etaname;
00195       etaname << "jet_eta_pmratio_" << i+1;
00196       histogramFactory().divide(histoPath(etaname.str()), *_h_eta_jet_plus[i], *_h_eta_jet_minus[i]);
00197       stringstream rapname;
00198       rapname << "jet_y_pmratio_" << i+1;
00199       histogramFactory().divide(histoPath(rapname.str()), *_h_rap_jet_plus[i], *_h_rap_jet_minus[i]);
00200     }
00201 
00202     for (int ibin = 0; ibin < _h_log10_R[m_njet]->size(); ++ibin) {
00203       IDataPoint* dp =_h_log10_R[m_njet]->point(ibin);
00204       dp->coordinate(1)->setValue(dp->coordinate(1)->value()*crossSection()/sumOfWeights());
00205     }
00206 
00207     // Scale the d{eta,R} histograms
00208     map<pair<size_t, size_t>, AIDA::IHistogram1D*>::iterator it;
00209     for (it = _h_deta_jets.begin(); it != _h_deta_jets.end(); ++it) {
00210       scale(it->second, crossSection()/sumOfWeights());
00211     }
00212     for (it = _h_dR_jets.begin(); it != _h_dR_jets.end(); ++it) {
00213       scale(it->second, crossSection()/sumOfWeights());
00214     }
00215 
00216     // Fill inclusive jet multi ratio
00217     int Nbins = _h_jet_multi_inclusive->axis().bins();
00218     std::vector<double> ratio(Nbins-1, 0.0);
00219     std::vector<double> err(Nbins-1, 0.0);
00220     for (int i = 0; i < Nbins-1; ++i) {
00221       if (_h_jet_multi_inclusive->binHeight(i) > 0.0 && _h_jet_multi_inclusive->binHeight(i+1) > 0.0) {
00222         ratio[i] = _h_jet_multi_inclusive->binHeight(i+1)/_h_jet_multi_inclusive->binHeight(i);
00223         double relerr_i = _h_jet_multi_inclusive->binError(i)/_h_jet_multi_inclusive->binHeight(i);
00224         double relerr_j = _h_jet_multi_inclusive->binError(i+1)/_h_jet_multi_inclusive->binHeight(i+1);
00225         err[i] = ratio[i] * (relerr_i + relerr_j);
00226       }
00227     }
00228     _h_jet_multi_ratio->setCoordinate(1, ratio, err);
00229 
00230     scale(_h_jet_multi_exclusive, crossSection()/sumOfWeights());
00231     scale(_h_jet_multi_inclusive, crossSection()/sumOfWeights());
00232   }
00233 
00234 
00235 }