MC_LEADINGJETS.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Tools/Logging.hh"
00004 #include "Rivet/RivetAIDA.hh"
00005 #include "Rivet/Projections/FinalState.hh"
00006 #include "Rivet/Projections/ChargedFinalState.hh"
00007 #include "Rivet/Projections/FastJets.hh"
00008 
00009 namespace Rivet {
00010 
00011 
00012   /* Underlying event in leading jet, extended to the LHC
00013    * @author Andy Buckley
00014    */
00015   class MC_LEADINGJETS : public Analysis {
00016   public:
00017  
00018     /// Constructor
00019     MC_LEADINGJETS()
00020       : Analysis("MC_LEADINGJETS")
00021     {    }
00022  
00023  
00024     /// @name Analysis methods
00025     //@{
00026  
00027     // Book histograms
00028     void init() {
00029       // Final state for the jet finding
00030       const FinalState fsj(-4.0, 4.0, 0.0*GeV);
00031       addProjection(fsj, "FSJ");
00032       addProjection(FastJets(fsj, FastJets::KT, 0.7), "Jets");
00033    
00034       // Charged final state for the distributions
00035       const ChargedFinalState cfs(-1.0, 1.0, 0.5*GeV);
00036       addProjection(cfs, "CFS");
00037 
00038       const double maxpt1 = 500.0;
00039       _hist_pnchg      = bookProfile1D("trans-nchg", 50, 0.0, maxpt1);
00040       _hist_pmaxnchg   = bookProfile1D("trans-maxnchg", 50, 0.0, maxpt1);
00041       _hist_pminnchg   = bookProfile1D("trans-minnchg", 50, 0.0, maxpt1);
00042       _hist_pcptsum    = bookProfile1D("trans-ptsum", 50, 0.0, maxpt1);
00043       _hist_pmaxcptsum = bookProfile1D("trans-maxptsum", 50, 0.0, maxpt1);
00044       _hist_pmincptsum = bookProfile1D("trans-minptsum", 50, 0.0, maxpt1);
00045       _hist_pcptave    = bookProfile1D("trans-ptavg", 50, 0.0, maxpt1);
00046     }
00047 
00048 
00049     // Do the analysis
00050     void analyze(const Event& e) {
00051 
00052       const FinalState& fsj = applyProjection<FinalState>(e, "FSJ");
00053       if (fsj.particles().empty()) {
00054         getLog() << Log::DEBUG << "Failed multiplicity cut" << endl;
00055         vetoEvent;
00056       }
00057 
00058       const FastJets& jetpro = applyProjection<FastJets>(e, "Jets");
00059       const Jets jets = jetpro.jetsByPt();
00060       getLog() << Log::DEBUG << "Jet multiplicity = " << jets.size() << endl;
00061 
00062       // Require the leading jet to be within |eta| < 2
00063       if (jets.size() < 1 || fabs(jets[0].momentum().pseudorapidity()) > 2) {
00064         getLog() << Log::DEBUG << "Failed jet cut" << endl;
00065         vetoEvent;
00066       }
00067 
00068       const double jetphi = jets[0].momentum().phi();
00069       const double jetpT  = jets[0].momentum().pT();
00070       getLog() << Log::DEBUG << "Leading jet: pT = " << jetpT/GeV << " GeV"
00071                << ", eta = " << jets[0].momentum().pseudorapidity()
00072                << ", phi = " << jetphi << endl;
00073 
00074       // Get the event weight
00075       const double weight = e.weight();
00076 
00077       // Get the final states to work with for filling the distributions
00078       const FinalState& cfs = applyProjection<ChargedFinalState>(e, "CFS");
00079 
00080       size_t   numOverall(0),     numToward(0),     numTrans1(0),     numTrans2(0),     numAway(0);
00081       double ptSumOverall(0.0), ptSumToward(0.0), ptSumTrans1(0.0), ptSumTrans2(0.0), ptSumAway(0.0);
00082       double ptMaxOverall(0.0), ptMaxToward(0.0), ptMaxTrans1(0.0), ptMaxTrans2(0.0), ptMaxAway(0.0);
00083 
00084       // Calculate all the charged stuff
00085       foreach (const Particle& p, cfs.particles()) {
00086         const double dPhi = deltaPhi(p.momentum().phi(), jetphi);
00087         const double pT = p.momentum().pT();
00088         const double phi = p.momentum().azimuthalAngle();
00089         const double rotatedphi = phi - jetphi;
00090 
00091         ptSumOverall += pT;
00092         ++numOverall;
00093         if (pT > ptMaxOverall) ptMaxOverall = pT;
00094 
00095         if (dPhi < PI/3.0) {
00096           ptSumToward += pT;
00097           ++numToward;
00098           if (pT > ptMaxToward) ptMaxToward = pT;
00099         }
00100         else if (dPhi < 2*PI/3.0) {
00101           if (rotatedphi <= PI) {
00102             ptSumTrans1 += pT;
00103             ++numTrans1;
00104             if (pT > ptMaxTrans1) ptMaxTrans1 = pT;
00105           } else {
00106             ptSumTrans2 += pT;
00107             ++numTrans2;
00108             if (pT > ptMaxTrans2) ptMaxTrans2 = pT;
00109           }
00110         }
00111         else {
00112           ptSumAway += pT;
00113           ++numAway;
00114           if (pT > ptMaxAway) ptMaxAway = pT;
00115         }
00116       }
00117    
00118    
00119       // Fill the histograms
00120       //_hist_tnchg->fill(jetpT/GeV, numToward/(4*PI/3), weight);
00121       _hist_pnchg->fill(jetpT/GeV, (numTrans1+numTrans2)/(4*PI/3), weight);
00122       _hist_pmaxnchg->fill(jetpT/GeV, (numTrans1>numTrans2 ? numTrans1 : numTrans2)/(2*PI/3), weight);
00123       _hist_pminnchg->fill(jetpT/GeV, (numTrans1<numTrans2 ? numTrans1 : numTrans2)/(2*PI/3), weight);
00124       //_hist_pdifnchg->fill(jetpT/GeV, abs(numTrans1-numTrans2)/(2*PI/3), weight);
00125       //_hist_anchg->fill(jetpT/GeV, numAway/(4*PI/3), weight);
00126    
00127       //_hist_tcptsum->fill(jetpT/GeV, ptSumToward/GeV/(4*PI/3), weight);
00128       _hist_pcptsum->fill(jetpT/GeV, (ptSumTrans1+ptSumTrans2)/GeV/(4*PI/3), weight);
00129       _hist_pmaxcptsum->fill(jetpT/GeV, (ptSumTrans1>ptSumTrans2 ? ptSumTrans1 : ptSumTrans2)/GeV/(2*PI/3), weight);
00130       _hist_pmincptsum->fill(jetpT/GeV, (ptSumTrans1<ptSumTrans2 ? ptSumTrans1 : ptSumTrans2)/GeV/(2*PI/3), weight);
00131       //_hist_pdifcptsum->fill(jetpT/GeV, fabs(ptSumTrans1-ptSumTrans2)/GeV/(2*PI/3), weight);
00132       //_hist_acptsum->fill(jetpT/GeV, ptSumAway/GeV/(4*PI/3), weight);
00133    
00134       //if (numToward > 0) {
00135       //  _hist_tcptave->fill(jetpT/GeV, ptSumToward/GeV/numToward, weight);
00136       //  _hist_tcptmax->fill(jetpT/GeV, ptMaxToward/GeV, weight);
00137       //}
00138       if ((numTrans1+numTrans2) > 0) {
00139         _hist_pcptave->fill(jetpT/GeV, (ptSumTrans1+ptSumTrans2)/GeV/(numTrans1+numTrans2), weight);
00140         //_hist_pcptmax->fill(jetpT/GeV, (ptMaxTrans1 > ptMaxTrans2 ? ptMaxTrans1 : ptMaxTrans2)/GeV, weight);
00141       }
00142       //if (numAway > 0) {
00143       //  _hist_acptave->fill(jetpT/GeV, ptSumAway/GeV/numAway, weight);
00144       //  _hist_acptmax->fill(jetpT/GeV, ptMaxAway/GeV, weight);
00145       //}
00146     }
00147  
00148  
00149     void finalize() {
00150       //
00151     }
00152  
00153  
00154   private:
00155  
00156     AIDA::IProfile1D *_hist_pnchg;
00157     AIDA::IProfile1D *_hist_pmaxnchg;
00158     AIDA::IProfile1D *_hist_pminnchg;
00159     AIDA::IProfile1D *_hist_pcptsum;
00160     AIDA::IProfile1D *_hist_pmaxcptsum;
00161     AIDA::IProfile1D *_hist_pmincptsum;
00162     AIDA::IProfile1D *_hist_pcptave;
00163  
00164   };
00165 
00166 
00167 
00168   // This global object acts as a hook for the plugin system
00169   AnalysisBuilder<MC_LEADINGJETS> plugin_MC_LEADINGJETS;
00170 
00171 }