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