MC_JetSplittings.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analyses/MC_JetSplittings.hh" 00003 #include "Rivet/Projections/FastJets.hh" 00004 00005 namespace Rivet { 00006 00007 00008 MC_JetSplittings::MC_JetSplittings(const string& name, 00009 size_t njet, 00010 const string& jetpro_name) 00011 : Analysis(name), m_njet(njet), m_jetpro_name(jetpro_name), 00012 _h_log10_d(njet), _h_log10_R(njet+1) 00013 { 00014 setNeedsCrossSection(true); // legitimate use, since a base class has no .info file! 00015 } 00016 00017 00018 00019 // Book histograms 00020 void MC_JetSplittings::init() { 00021 for (size_t i = 0; i < m_njet; ++i) { 00022 string dname = "log10_d_" + to_str(i) + to_str(i+1); 00023 _h_log10_d[i] = bookHisto1D(dname, 100, 0.2, log10(0.5*sqrtS())); 00024 string Rname = "log10_R_" + to_str(i); 00025 _h_log10_R[i] = bookScatter2D(Rname, 50, 0.2, log10(0.5*sqrtS())); 00026 } 00027 string Rname = "log10_R_" + to_str(m_njet); 00028 _h_log10_R[m_njet] = bookScatter2D(Rname, 50, 0.2, log10(0.5*sqrtS())); 00029 } 00030 00031 00032 00033 // Do the analysis 00034 void MC_JetSplittings::analyze(const Event & e) { 00035 const double weight = e.weight(); 00036 00037 const FastJets& jetpro = applyProjection<FastJets>(e, m_jetpro_name); 00038 00039 // Jet resolutions and integrated jet rates 00040 const fastjet::ClusterSequence* seq = jetpro.clusterSeq(); 00041 if (seq != NULL) { 00042 double previous_dij = 10.0; 00043 for (size_t i = 0; i < min(m_njet,(size_t)seq->n_particles()); ++i) { 00044 // Jet resolution i -> j 00045 double d_ij = log10(sqrt(seq->exclusive_dmerge_max(i))); 00046 00047 // Fill differential jet resolution 00048 _h_log10_d[i]->fill(d_ij, weight); 00049 00050 // Fill integrated jet resolution 00051 for (size_t ibin = 0; ibin < _h_log10_R[i]->numPoints(); ++ibin) { 00052 Point2D& dp = _h_log10_R[i]->point(ibin); 00053 if (dp.x() > d_ij && dp.x() < previous_dij) { 00054 dp.setY(dp.y() + weight); 00055 } 00056 } 00057 previous_dij = d_ij; 00058 } 00059 // One remaining integrated jet resolution 00060 for (size_t ibin = 0; ibin<_h_log10_R[m_njet]->numPoints(); ++ibin) { 00061 Point2D & dp = _h_log10_R[m_njet]->point(ibin); 00062 if (dp.x() < previous_dij) { 00063 dp.setY(dp.y() + weight); 00064 } 00065 } 00066 } 00067 00068 } 00069 00070 00071 // Finalize 00072 void MC_JetSplittings::finalize() { 00073 for (size_t i = 0; i < m_njet; ++i) { 00074 scale(_h_log10_d[i], crossSection()/sumOfWeights()); 00075 for (size_t ibin = 0; ibin<_h_log10_R[i]->numPoints(); ++ibin) { 00076 Point2D & dp = _h_log10_R[i]->point(ibin); 00077 dp.setY(dp.y()*crossSection()/sumOfWeights()); 00078 } 00079 } 00080 00081 for (size_t ibin = 0; ibin < _h_log10_R[m_njet]->numPoints(); ++ibin) { 00082 Point2D & dp =_h_log10_R[m_njet]->point(ibin); 00083 dp.setY(dp.y()*crossSection()/sumOfWeights()); 00084 } 00085 } 00086 00087 00088 } Generated on Fri Oct 25 2013 12:41:46 for The Rivet MC analysis system by ![]() |