D0_2008_S6879055.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/RivetAIDA.hh"
00004 #include "Rivet/Tools/Logging.hh"
00005 #include "Rivet/Projections/FinalState.hh"
00006 #include "Rivet/Projections/LeadingParticlesFinalState.hh"
00007 #include "Rivet/Projections/InvMassFinalState.hh"
00008 #include "Rivet/Projections/VetoedFinalState.hh"
00009 #include "Rivet/Projections/FastJets.hh"
00010 
00011 namespace Rivet {
00012 
00013 
00014   /// @brief Measurement of the ratio sigma(Z/gamma* + n jets)/sigma(Z/gamma*)
00015   class D0_2008_S6879055 : public Analysis {
00016   public:
00017  
00018     /// Default constructor.
00019     D0_2008_S6879055() : Analysis("D0_2008_S6879055")
00020     {
00021       setBeams(PROTON, ANTIPROTON);
00022     }
00023 
00024 
00025     /// @name Analysis methods
00026     //@{
00027  
00028     // Book histograms
00029     void init() {
00030       // Basic final state
00031       FinalState fs(-5.0, 5.0);
00032       addProjection(fs, "FS");
00033    
00034       // Leading electrons in tracking acceptance
00035       LeadingParticlesFinalState lpfs(FinalState(-1.1, 1.1, 25*GeV));
00036       lpfs.addParticleId(ELECTRON).addParticleId(POSITRON);
00037       addProjection(lpfs, "LeadingElectronsFS");
00038    
00039       // Invariant mass selection around Z pole
00040       InvMassFinalState electronsFromZ(lpfs, make_pair(ELECTRON, POSITRON), 75*GeV, 105*GeV);
00041       addProjection(electronsFromZ,"ElectronsFromZ");
00042    
00043       // Vetoed FS for jets
00044       VetoedFinalState vfs(fs);
00045       // Add particle/antiparticle vetoing
00046       vfs.vetoNeutrinos();
00047       // Veto the electrons from Z decay
00048       vfs.addVetoOnThisFinalState(electronsFromZ);
00049       addProjection(vfs, "JetFS");
00050    
00051       // Jet finder
00052       FastJets jets(vfs, FastJets::D0ILCONE, 0.5);
00053       addProjection(jets, "Jets");
00054    
00055       _crossSectionRatio = bookHistogram1D(1, 1, 1);
00056       _pTjet1 = bookHistogram1D(2, 1, 1);
00057       _pTjet2 = bookHistogram1D(3, 1, 1);
00058       _pTjet3 = bookHistogram1D(4, 1, 1);
00059     }
00060  
00061  
00062  
00063     /// Do the analysis
00064     void analyze(const Event& event) {   
00065       // Skip if the event is empty
00066       const FinalState& fs = applyProjection<FinalState>(event, "FS");
00067       if (fs.empty()) vetoEvent;
00068       const double weight = event.weight();
00069    
00070       // Find the Z candidates
00071       const InvMassFinalState& invmassfs = applyProjection<InvMassFinalState>(event, "ElectronsFromZ");
00072       // If there is no Z candidate in the FinalState, skip the event
00073       if (invmassfs.particles().size() != 2) {
00074         getLog() << Log::DEBUG << "No Z candidate found" << endl;
00075         vetoEvent;
00076       }
00077    
00078       // Now build the list of jets on a FS without the electrons from Z
00079       // Additional cuts on jets: |eta| < 2.5 and dR(j,leading electron) > 0.4
00080       const JetAlg& jetpro = applyProjection<JetAlg>(event, "Jets");
00081       const Jets jets = jetpro.jetsByPt(20.0*GeV);
00082       vector<FourMomentum> finaljet_list;
00083       foreach (const Jet& j, jets) {
00084         const double jeta = j.momentum().eta();
00085         const double jphi = j.momentum().phi();
00086         if (fabs(jeta) > 2.5) continue;
00087      
00088         FourMomentum e0 = invmassfs.particles()[0].momentum();
00089         FourMomentum e1 = invmassfs.particles()[1].momentum();
00090         const double e0eta = e0.pseudorapidity();
00091         const double e0phi = e0.azimuthalAngle();
00092         if (deltaR(e0eta, e0phi, jeta, jphi) < 0.4) continue;
00093      
00094         const double e1eta = e1.pseudorapidity();
00095         const double e1phi = e1.azimuthalAngle();
00096         if (deltaR(e1eta, e1phi, jeta, jphi) < 0.4) continue;
00097      
00098         // If we pass all cuts...
00099         finaljet_list.push_back(j.momentum());
00100       }
00101       getLog() << Log::DEBUG << "Num jets passing = " << finaljet_list.size() << endl;
00102    
00103       // For normalisation of crossSection data (includes events with no jets passing cuts)
00104       _crossSectionRatio->fill(0, weight);
00105    
00106       // Fill jet pT and multiplicities
00107       if (finaljet_list.size() >= 1) {
00108         _crossSectionRatio->fill(1, weight);
00109         _pTjet1->fill(finaljet_list[0].pT(), weight);
00110       }
00111       if (finaljet_list.size() >= 2) {
00112         _crossSectionRatio->fill(2, weight);
00113         _pTjet2->fill(finaljet_list[1].pT(), weight);
00114       }
00115       if (finaljet_list.size() >= 3) {
00116         _crossSectionRatio->fill(3, weight);
00117         _pTjet3->fill(finaljet_list[2].pT(), weight);
00118       }
00119       if (finaljet_list.size() >= 4) {
00120         _crossSectionRatio->fill(4, weight);
00121       }
00122     }
00123  
00124  
00125  
00126     /// Finalize
00127     void finalize() {
00128       // Now divide by the inclusive result
00129       _crossSectionRatio->scale(1.0/_crossSectionRatio->binHeight(0));
00130    
00131       // Normalise jet pTs to integrals of data
00132       // NB. There is no other way to do this, because these quantities are not
00133       // detector-corrected
00134       normalize(_pTjet1, 10439.0); // fixed norm OK
00135       normalize(_pTjet2, 1461.5); // fixed norm OK
00136       normalize(_pTjet3, 217.0); // fixed norm OK
00137     }
00138  
00139     //@}
00140 
00141 
00142   private:
00143 
00144     /// @name Histograms
00145     //@{
00146     AIDA::IHistogram1D * _crossSectionRatio;
00147     AIDA::IHistogram1D * _pTjet1;
00148     AIDA::IHistogram1D * _pTjet2;
00149     AIDA::IHistogram1D * _pTjet3;
00150     //@}
00151 
00152   };
00153 
00154  
00155  
00156   // This global object acts as a hook for the plugin system
00157   AnalysisBuilder<D0_2008_S6879055> plugin_D0_2008_S6879055;
00158 
00159 }