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