rivet is hosted by Hepforge, IPPP Durham
ATLAS_2010_CONF_2010_049.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/RivetYODA.hh"
00004 #include "Rivet/Projections/ChargedFinalState.hh"
00005 #include "Rivet/Projections/FastJets.hh"
00006 #include "Rivet/Tools/Logging.hh"
00007 
00008 namespace Rivet {
00009 
00010 
00011   class ATLAS_2010_CONF_2010_049 : public Analysis {
00012   public:
00013 
00014     ATLAS_2010_CONF_2010_049()
00015       : Analysis("ATLAS_2010_CONF_2010_049")
00016     {    }
00017 
00018 
00019     void init() {
00020       ChargedFinalState cfs(-1.5, 1.5, 0.5*GeV);
00021       addProjection(cfs, "CFS");
00022 
00023       FastJets jetsproj6(cfs, FastJets::ANTIKT, 0.6);
00024       addProjection(jetsproj6, "Jets6");
00025 
00026       FastJets jetsproj4(cfs, FastJets::ANTIKT, 0.4);
00027       addProjection(jetsproj4, "Jets4");
00028 
00029       for (size_t i=0 ; i<2 ; i++) {
00030         _h_xsec[i]       = bookHisto1D(1+i, 1, 1);
00031         _h_frag_04_06[i] = bookHisto1D(3+i, 1, 1);
00032         _h_frag_06_10[i] = bookHisto1D(3+i, 2, 1);
00033         _h_frag_10_15[i] = bookHisto1D(3+i, 3, 1);
00034         _h_frag_15_24[i] = bookHisto1D(3+i, 4, 1);
00035         _njets_04_06[i] = 0.0;
00036         _njets_06_10[i] = 0.0;
00037         _njets_10_15[i] = 0.0;
00038         _njets_15_24[i] = 0.0;
00039       }
00040     }
00041 
00042 
00043     void analyze(const Event& event) {
00044       const double weight = event.weight();
00045 
00046       const FastJets & jetsproj6 = applyProjection<FastJets>(event, "Jets6");
00047       const FastJets & jetsproj4 = applyProjection<FastJets>(event, "Jets4");
00048       Jets alljets[2];
00049       alljets[0] = jetsproj6.jetsByPt(4.0*GeV);
00050       alljets[1] = jetsproj4.jetsByPt(4.0*GeV);
00051 
00052       for (size_t i=0 ; i<2 ; i++) {
00053         Jets jets;
00054 
00055         // First we want to make sure that we only use jets within |eta|<0.57
00056         foreach (const Jet& jet, alljets[i]) {
00057           if (fabs(jet.momentum().eta())<0.57) {
00058             jets.push_back(jet);
00059           }
00060         }
00061         foreach (const Jet& jet, jets) {
00062           const double pTjet = jet.momentum().pT();
00063           const double pjet = jet.momentum().p().mod();
00064           _h_xsec[i]->fill(pTjet, weight);
00065           if (pTjet > 24*GeV) continue;
00066           foreach (const Particle& p, jet.particles()) {
00067             double z=p.momentum().p().mod()/pjet;
00068             if (z>0.9999) z=0.9999;   // Make sure that z=1 doesn't go into overflow
00069             if (pTjet > 15*GeV) {
00070               _h_frag_15_24[i]->fill(z, weight);
00071             }
00072             else if (pTjet > 10*GeV) {
00073               _h_frag_10_15[i]->fill(z, weight);
00074             }
00075             else if (pTjet > 6*GeV) {
00076               _h_frag_06_10[i]->fill(z, weight);
00077             }
00078             else {
00079               _h_frag_04_06[i]->fill(z, weight);
00080             }
00081           }
00082           if (pTjet > 15*GeV) {
00083             _njets_15_24[i] += weight;
00084           }
00085           else if (pTjet > 10*GeV) {
00086             _njets_10_15[i] += weight;
00087           }
00088           else if (pTjet > 6*GeV) {
00089             _njets_06_10[i] += weight;
00090           }
00091           else {
00092             _njets_04_06[i] += weight;
00093           }
00094         }
00095       }
00096     }
00097 
00098     void finalize() {
00099       for (size_t i=0 ; i<2 ; i++) {
00100         // deta = 2*0.57
00101         scale(_h_xsec[i], crossSection()/microbarn/sumOfWeights()/(2*0.57));
00102         scale(_h_frag_04_06[i], 1./_njets_04_06[i]);
00103         scale(_h_frag_06_10[i], 1./_njets_06_10[i]);
00104         scale(_h_frag_10_15[i], 1./_njets_10_15[i]);
00105         scale(_h_frag_15_24[i], 1./_njets_15_24[i]);
00106       }
00107     }
00108 
00109 
00110   private:
00111 
00112     Histo1DPtr _h_xsec[2];
00113     Histo1DPtr _h_frag_04_06[2];
00114     Histo1DPtr _h_frag_06_10[2];
00115     Histo1DPtr _h_frag_10_15[2];
00116     Histo1DPtr _h_frag_15_24[2];
00117     double _njets_04_06[2];
00118     double _njets_06_10[2];
00119     double _njets_10_15[2];
00120     double _njets_15_24[2];
00121   };
00122 
00123 
00124 
00125   // The hook for the plugin system
00126   DECLARE_RIVET_PLUGIN(ATLAS_2010_CONF_2010_049);
00127 
00128 }