rivet is hosted by Hepforge, IPPP Durham
CDF_2001_S4751469.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 
00003 // Field & Stuart underlying event analysis at CDF.
00004 // Phys.Rev.D65:092002,2002 - no arXiv code.
00005 // FNAL-PUB 01/211-E
00006 
00007 #include "Rivet/Analysis.hh"
00008 #include "Rivet/Projections/ChargedFinalState.hh"
00009 #include "Rivet/Projections/ConstLossyFinalState.hh"
00010 #include "Rivet/Projections/FastJets.hh"
00011 #include "Rivet/Projections/TriggerCDFRun0Run1.hh"
00012 
00013 namespace Rivet {
00014 
00015 
00016   /// @brief "Field-Stuart" CDF Run I track-jet underlying event analysis
00017   ///
00018   /// @author Andy Buckley
00019   ///
00020   /// The "original" underlying event analysis, using a non-standard track-jet algorithm.
00021   ///
00022   /// @par Run conditions
00023   ///
00024   /// @arg \f$ \sqrt{s} = \f$ 1800 GeV
00025   /// @arg Run with generic QCD events.
00026   /// @arg Several \f$ p_\perp^\text{min} \f$ cutoffs are probably required to fill the profile histograms:
00027   ///   @arg \f$ p_\perp^\text{min} = \f$ 0 (min bias), 10, 20 GeV
00028   class CDF_2001_S4751469 : public Analysis {
00029   public:
00030 
00031     /// Constructor: cuts on final state are \f$ -1 < \eta < 1 \f$
00032     /// and \f$ p_T > 0.5 \f$ GeV.
00033     CDF_2001_S4751469()
00034       : Analysis("CDF_2001_S4751469"),
00035         _totalNumTrans2(0), _totalNumTrans5(0), _totalNumTrans30(0),
00036         _sumWeightsPtLead2(0),_sumWeightsPtLead5(0), _sumWeightsPtLead30(0)
00037     {
00038     }
00039 
00040 
00041     /// @name Analysis methods
00042     //@{
00043 
00044     // Book histograms
00045     void init() {
00046       addProjection(TriggerCDFRun0Run1(), "Trigger");
00047       // Randomly discard 8% of charged particles as a kind of hacky detector correction.
00048       const ChargedFinalState cfs(-1.0, 1.0, 0.5*GeV);
00049       const ConstLossyFinalState lfs(cfs, 0.08);
00050 
00051       addProjection(lfs, "FS");
00052       addProjection(FastJets(lfs, FastJets::TRACKJET, 0.7), "TrackJet");
00053 
00054       _numvsDeltaPhi2 =  bookProfile1D(1, 1, 1);
00055       _numvsDeltaPhi5 =  bookProfile1D(1, 1, 2);
00056       _numvsDeltaPhi30 = bookProfile1D(1, 1, 3);
00057       _pTvsDeltaPhi2 =   bookProfile1D(2, 1, 1);
00058       _pTvsDeltaPhi5 =   bookProfile1D(2, 1, 2);
00059       _pTvsDeltaPhi30 =  bookProfile1D(2, 1, 3);
00060 
00061       _numTowardMB = bookProfile1D(3, 1, 1);
00062       _numTransMB = bookProfile1D(3, 1, 2);
00063       _numAwayMB = bookProfile1D(3, 1, 3);
00064       _numTowardJ20 = bookProfile1D(4, 1, 1);
00065       _numTransJ20 = bookProfile1D(4, 1, 2);
00066       _numAwayJ20 = bookProfile1D(4, 1, 3);
00067 
00068       _ptsumTowardMB = bookProfile1D(5, 1, 1);
00069       _ptsumTransMB = bookProfile1D(5, 1, 2);
00070       _ptsumAwayMB = bookProfile1D(5, 1, 3);
00071       _ptsumTowardJ20 = bookProfile1D(6, 1, 1);
00072       _ptsumTransJ20 = bookProfile1D(6, 1, 2);
00073       _ptsumAwayJ20 = bookProfile1D(6, 1, 3);
00074 
00075       _ptTrans2 = bookHisto1D(7, 1, 1);
00076       _ptTrans5 = bookHisto1D(7, 1, 2);
00077       _ptTrans30 = bookHisto1D(7, 1, 3);
00078     }
00079 
00080 
00081     /// Do the analysis
00082     void analyze(const Event& event) {
00083       // Trigger
00084       const bool trigger = applyProjection<TriggerCDFRun0Run1>(event, "Trigger").minBiasDecision();
00085       if (!trigger) vetoEvent;
00086 
00087       // Analyse, with pT > 0.5 GeV AND |eta| < 1
00088       const JetAlg& tj = applyProjection<JetAlg>(event, "TrackJet");
00089 
00090       // Final state (lossy) charged particles
00091       const FinalState& fs = applyProjection<FinalState>(event, "FS");
00092 
00093       // Get jets, sorted by pT
00094       const Jets jets = tj.jetsByPt();
00095       if (jets.empty()) {
00096         vetoEvent;
00097       }
00098 
00099       Jet leadingJet = jets.front();
00100       const double phiLead = leadingJet.phi();
00101       const double ptLead = leadingJet.pT();
00102 
00103       // Cut on highest pT jet: combined 0.5 GeV < pT(lead) < 50 GeV
00104       if (ptLead/GeV < 0.5) vetoEvent;
00105       if (ptLead/GeV > 50.0) vetoEvent;
00106 
00107       // Get the event weight
00108       const double weight = event.weight();
00109 
00110       // Count sum of all event weights in three pT_lead regions
00111       if (ptLead/GeV > 2.0) {
00112         _sumWeightsPtLead2 += weight;
00113       }
00114       if (ptLead/GeV > 5.0) {
00115         _sumWeightsPtLead5 += weight;
00116       }
00117       if (ptLead/GeV > 30.0) {
00118         _sumWeightsPtLead30 += weight;
00119       }
00120 
00121       // Run over tracks
00122       double ptSumToward(0.0), ptSumAway(0.0), ptSumTrans(0.0);
00123       size_t numToward(0), numTrans(0), numAway(0);
00124 
00125       // Temporary histos that bin N and pT in dphi
00126       Profile1D htmp_num_dphi_2(refData(1, 1, 1)), htmp_num_dphi_5(refData(1, 1, 2)), htmp_num_dphi_30(refData(1, 1, 3));
00127       Profile1D htmp_pt_dphi_2(refData(2, 1, 1)), htmp_pt_dphi_5(refData(2, 1, 2)), htmp_pt_dphi_30(refData(2, 1, 3));
00128 
00129       foreach (const Particle& p, fs.particles()) {
00130         // Calculate DeltaPhi(p,leadingJet)
00131         const double dPhi = deltaPhi(p.momentum().phi(), phiLead);
00132         const double pT = p.pT();
00133 
00134         if (dPhi < PI/3.0) {
00135           ptSumToward += pT;
00136           ++numToward;
00137         }
00138         else if (dPhi < 2*PI/3.0) {
00139           ptSumTrans += pT;
00140           ++numTrans;
00141           // Fill transverse pT distributions
00142           if (ptLead/GeV > 2.0) {
00143             _ptTrans2->fill(pT/GeV, weight);
00144             _totalNumTrans2 += weight;
00145           }
00146           if (ptLead/GeV > 5.0) {
00147             _ptTrans5->fill(pT/GeV, weight);
00148             _totalNumTrans5 += weight;
00149           }
00150           if (ptLead/GeV > 30.0) {
00151             _ptTrans30->fill(pT/GeV, weight);
00152             _totalNumTrans30 += weight;
00153           }
00154         }
00155         else {
00156           ptSumAway += pT;
00157           ++numAway;
00158         }
00159 
00160         // Fill tmp histos to bin event's track Nch & pT in dphi
00161         const double dPhideg = 180*dPhi/PI;
00162         if (ptLead/GeV > 2.0) {
00163           htmp_num_dphi_2.fill(dPhideg, 1);
00164           htmp_pt_dphi_2.fill (dPhideg, pT/GeV);
00165         }
00166         if (ptLead/GeV > 5.0) {
00167           htmp_num_dphi_5.fill(dPhideg, 1);
00168           htmp_pt_dphi_5.fill (dPhideg, pT/GeV);
00169         }
00170         if (ptLead/GeV > 30.0) {
00171           htmp_num_dphi_30.fill(dPhideg, 1);
00172           htmp_pt_dphi_30.fill (dPhideg, pT/GeV);
00173         }
00174       }
00175 
00176       // Update the "proper" dphi profile histograms
00177       for (int i = 0; i < 50; i++) { //< @todo Should really explicitly iterate over nbins for each temp histo
00178         if (ptLead/GeV > 2.0) {
00179           const double x2 = htmp_pt_dphi_2.bin(i).midpoint();
00180           const double num2 = (htmp_num_dphi_2.bin(i).numEntries() > 0) ? htmp_num_dphi_2.bin(i).mean() : 0.0;
00181           const double pt2 = (htmp_num_dphi_2.bin(i).numEntries() > 0) ? htmp_pt_dphi_2.bin(i).mean() : 0.0;
00182           _numvsDeltaPhi2->fill(x2, num2, weight);
00183           _pTvsDeltaPhi2->fill(x2, pt2, weight);
00184         }
00185         if (ptLead/GeV > 5.0) {
00186           const double x5 = htmp_pt_dphi_5.bin(i).midpoint();
00187           const double num5 = (htmp_num_dphi_5.bin(i).numEntries() > 0) ? htmp_num_dphi_5.bin(i).mean() : 0.0;
00188           const double pt5 = (htmp_num_dphi_5.bin(i).numEntries() > 0) ? htmp_pt_dphi_5.bin(i).mean() : 0.0;
00189           _numvsDeltaPhi5->fill(x5, num5, weight);
00190           _pTvsDeltaPhi5->fill(x5, pt5, weight);
00191         }
00192         if (ptLead/GeV > 30.0) {
00193           const double x30 = htmp_pt_dphi_30.bin(i).midpoint();
00194           const double num30 = (htmp_num_dphi_30.bin(i).numEntries() > 0) ? htmp_num_dphi_30.bin(i).mean() : 0.0;
00195           const double pt30 = (htmp_num_dphi_30.bin(i).numEntries() > 0) ? htmp_pt_dphi_30.bin(i).mean() : 0.0;
00196           _numvsDeltaPhi30->fill(x30, num30, weight);
00197           _pTvsDeltaPhi30->fill(x30, pt30, weight);
00198         }
00199       }
00200 
00201       // Log some event details about pT
00202       MSG_DEBUG("pT [lead; twd, away, trans] = [" << ptLead << "; "
00203                 << ptSumToward << ", " << ptSumAway << ", " << ptSumTrans << "]");
00204 
00205       // Update the pT profile histograms
00206       _ptsumTowardMB->fill(ptLead/GeV, ptSumToward/GeV, weight);
00207       _ptsumTowardJ20->fill(ptLead/GeV, ptSumToward/GeV, weight);
00208 
00209       _ptsumTransMB->fill(ptLead/GeV, ptSumTrans/GeV, weight);
00210       _ptsumTransJ20->fill(ptLead/GeV, ptSumTrans/GeV, weight);
00211 
00212       _ptsumAwayMB->fill(ptLead/GeV, ptSumAway/GeV, weight);
00213       _ptsumAwayJ20->fill(ptLead/GeV, ptSumAway/GeV, weight);
00214 
00215       // Log some event details about Nch
00216       MSG_DEBUG("N [twd, away, trans] = [" << ptLead << "; "
00217                 << numToward << ", " << numTrans << ", " << numAway << "]");
00218 
00219       // Update the N_track profile histograms
00220       _numTowardMB->fill(ptLead/GeV, numToward, weight);
00221       _numTowardJ20->fill(ptLead/GeV, numToward, weight);
00222 
00223       _numTransMB->fill(ptLead/GeV, numTrans, weight);
00224       _numTransJ20->fill(ptLead/GeV, numTrans, weight);
00225 
00226       _numAwayMB->fill(ptLead/GeV, numAway, weight);
00227       _numAwayJ20->fill(ptLead/GeV, numAway, weight);
00228     }
00229 
00230 
00231     /// Normalize histos
00232     void finalize() {
00233       normalize(_ptTrans2, _totalNumTrans2 / _sumWeightsPtLead2);
00234       normalize(_ptTrans5, _totalNumTrans5 / _sumWeightsPtLead5);
00235       normalize(_ptTrans30, _totalNumTrans30 / _sumWeightsPtLead30);
00236     }
00237 
00238     //@}
00239 
00240 
00241   private:
00242 
00243     /// Sum total number of charged particles in the trans region, in 3 \f$ p_\perp^\text{lead} \f$ bins.
00244     double _totalNumTrans2, _totalNumTrans5, _totalNumTrans30;
00245 
00246     /// Sum the total number of events in 3 \f$ p_\perp^\text{lead} \f$ bins.
00247     double _sumWeightsPtLead2,_sumWeightsPtLead5, _sumWeightsPtLead30;
00248 
00249 
00250     /// @name Histogram collections
00251     //@{
00252 
00253     // The sumpt vs. dphi and Nch vs. dphi histos
00254     Profile1DPtr _numvsDeltaPhi2, _numvsDeltaPhi5, _numvsDeltaPhi30;
00255     Profile1DPtr _pTvsDeltaPhi2, _pTvsDeltaPhi5, _pTvsDeltaPhi30;
00256 
00257     /// Profile histograms, binned in the \f$ p_T \f$ of the leading jet, for
00258     /// the \f$ p_T \f$ sum in the toward, transverse and away regions.
00259     Profile1DPtr _ptsumTowardMB,  _ptsumTransMB,  _ptsumAwayMB;
00260     Profile1DPtr _ptsumTowardJ20, _ptsumTransJ20, _ptsumAwayJ20;
00261 
00262     /// Profile histograms, binned in the \f$ p_T \f$ of the leading jet, for
00263     /// the number of charged particles per jet in the toward, transverse and
00264     /// away regions.
00265     Profile1DPtr _numTowardMB,  _numTransMB,  _numAwayMB;
00266     Profile1DPtr _numTowardJ20, _numTransJ20, _numAwayJ20;
00267 
00268     /// Histogram of \f$ p_T \f$ distribution for 3 different \f$ p_{T1} \f$ IR cutoffs.
00269     Histo1DPtr _ptTrans2, _ptTrans5, _ptTrans30;
00270     //@}
00271 
00272   };
00273 
00274 
00275 
00276   // The hook for the plugin system
00277   DECLARE_RIVET_PLUGIN(CDF_2001_S4751469);
00278 
00279 }