CDF_2008_S8095620.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/Tools/ParticleIdUtils.hh"
00006 #include "Rivet/Projections/FastJets.hh"
00007 #include "Rivet/Projections/FinalState.hh"
00008 #include "Rivet/Projections/VetoedFinalState.hh"
00009 #include "Rivet/Projections/InvMassFinalState.hh"
00010 
00011 namespace Rivet {
00012 
00013 
00014   /// @brief CDF Run II Z + b-jet cross-section measurement
00015   class CDF_2008_S8095620 : public Analysis {
00016   public:
00017     
00018     /// Constructor.
00019     /// jet cuts: |eta| <= 1.5
00020     CDF_2008_S8095620()
00021       : Analysis("CDF_2008_S8095620"),
00022         _Rjet(0.7), _JetPtCut(20.), _JetEtaCut(1.5), _Lep1PtCut(18.), _Lep2PtCut(10.), _LepEtaCut(3.2),
00023         _sumWeightSelected(0.0)
00024     {
00025       setBeams(PROTON, ANTIPROTON);
00026     }
00027  
00028 
00029     /// @name Analysis methods
00030     //@{
00031  
00032     void init() {
00033       // Set up projections
00034       const FinalState fs(-3.2, 3.2);
00035       addProjection(fs, "FS");
00036       // Create a final state with any e+e- or mu+mu- pair with
00037       // invariant mass 76 -> 106 GeV and ET > 18 (Z decay products)
00038       vector<pair<PdgId,PdgId> > vids;
00039       vids.push_back(make_pair(ELECTRON, POSITRON));
00040       vids.push_back(make_pair(MUON, ANTIMUON));
00041       FinalState fs2(-3.2, 3.2);
00042       InvMassFinalState invfs(fs2, vids, 76*GeV, 106*GeV);
00043       addProjection(invfs, "INVFS");
00044       // Make a final state without the Z decay products for jet clustering
00045       VetoedFinalState vfs(fs);
00046       vfs.addVetoOnThisFinalState(invfs);
00047       addProjection(vfs, "VFS");
00048       addProjection(FastJets(vfs, FastJets::CDFMIDPOINT, 0.7), "Jets");
00049 
00050       // Book histograms
00051       _dStot    = bookHistogram1D(1, 1, 1);
00052       _dSdET    = bookHistogram1D(2, 1, 1);
00053       _dSdETA   = bookHistogram1D(3, 1, 1);
00054       _dSdZpT   = bookHistogram1D(4, 1, 1);
00055       _dSdNJet  = bookHistogram1D(5, 1, 1);
00056       _dSdNbJet = bookHistogram1D(6, 1, 1);
00057      }
00058  
00059 
00060     // Do the analysis
00061     void analyze(const Event& event) {
00062       // Check we have an l+l- pair that passes the kinematic cuts
00063       // Get the Z decay products (mu+mu- or e+e- pair)
00064       const InvMassFinalState& invMassFinalState = applyProjection<InvMassFinalState>(event, "INVFS");
00065       const ParticleVector&  ZDecayProducts =  invMassFinalState.particles();
00066    
00067       // make sure we have 2 Z decay products (mumu or ee)
00068       if (ZDecayProducts.size() < 2) vetoEvent;
00069       //new cuts
00070       double Lep1Pt = ZDecayProducts[0].momentum().perp();
00071       double Lep2Pt = ZDecayProducts[1].momentum().perp();
00072       double Lep1Eta = fabs(ZDecayProducts[0].momentum().rapidity());
00073       double Lep2Eta = fabs(ZDecayProducts[1].momentum().rapidity());
00074 
00075       if (Lep1Eta > _LepEtaCut || Lep2Eta > _LepEtaCut) vetoEvent;
00076 
00077       if (abs(ZDecayProducts[0].pdgId())==13 && 
00078           ((Lep1Eta > 1.5 || Lep2Eta > 1.5) || (Lep1Eta > 1.0 && Lep2Eta > 1.0))) {
00079         vetoEvent;
00080       }
00081       
00082       if (Lep1Pt > Lep2Pt) {
00083         if (Lep1Pt < _Lep1PtCut || Lep2Pt < _Lep2PtCut) vetoEvent;
00084       }
00085       else {
00086         if (Lep1Pt < _Lep2PtCut || Lep2Pt < _Lep1PtCut) vetoEvent;
00087       }
00088       
00089       _sumWeightSelected += event.weight();
00090       // @todo: write out a warning if there are more than two decay products
00091       FourMomentum Zmom = ZDecayProducts[0].momentum() +  ZDecayProducts[1].momentum();
00092    
00093       // Put all b-quarks in a vector
00094       ParticleVector bquarks;
00095       foreach (const GenParticle* p, particles(event.genEvent())) {
00096         if (fabs(p->pdg_id()) == BQUARK) {
00097           bquarks += Particle(*p);
00098         }
00099       }
00100    
00101       // Get jets
00102       const FastJets& jetpro = applyProjection<FastJets>(event, "Jets");
00103       getLog() << Log::DEBUG << "Jet multiplicity before any pT cut = " << jetpro.size() << endl;
00104    
00105       const PseudoJets& jets = jetpro.pseudoJetsByPt();
00106       getLog() << Log::DEBUG << "jetlist size = " << jets.size() << endl;
00107    
00108       int numBJet = 0;
00109       int numJet  = 0;
00110       // for each b-jet plot the ET and the eta of the jet, normalise to the total cross section at the end
00111       // for each event plot N jet and pT(Z), normalise to the total cross section at the end
00112       for (PseudoJets::const_iterator jt = jets.begin(); jt != jets.end(); ++jt) {
00113         // select jets that pass the kinematic cuts
00114         if (jt->perp() > _JetPtCut && fabs(jt->rapidity()) <= _JetEtaCut) {
00115           numJet++;
00116           // does the jet contain a b-quark?
00117           bool bjet = false;
00118           foreach (const Particle& bquark,  bquarks) {
00119             if (deltaR(jt->rapidity(), jt->phi(), bquark.momentum().rapidity(),bquark.momentum().azimuthalAngle()) <= _Rjet) {
00120               bjet = true;
00121               break;
00122             }
00123           } // end loop around b-jets
00124           if (bjet) {
00125             numBJet++;
00126             _dSdET->fill(jt->perp(),event.weight());
00127             _dSdETA->fill(fabs(jt->rapidity()),event.weight());
00128           }
00129         }
00130       } // end loop around jets
00131    
00132       // wasn't asking for b-jets before!!!!
00133       if(numJet > 0 && numBJet > 0) _dSdNJet->fill(numJet,event.weight());
00134       if(numBJet > 0) {
00135         _dStot->fill(1960.0,event.weight());
00136         _dSdNbJet->fill(numBJet,event.weight());
00137         _dSdZpT->fill(Zmom.pT(),event.weight());
00138       }
00139     }
00140  
00141 
00142 
00143     // Finalize
00144     void finalize() {
00145       // normalise histograms
00146       // scale by 1 / the sum-of-weights of events that pass the Z cuts
00147       // since the cross sections are normalized to the inclusive
00148       // Z cross sections.
00149       double Scale = 1.0;
00150       if (_sumWeightSelected != 0.0) Scale = 1.0/_sumWeightSelected;
00151       _dStot->scale(Scale);
00152       _dSdET->scale(Scale);
00153       _dSdETA->scale(Scale);
00154       _dSdNJet->scale(Scale);
00155       _dSdNbJet->scale(Scale);
00156       _dSdZpT->scale(Scale);
00157     }
00158 
00159     //@}
00160 
00161 
00162   private:
00163 
00164     double _Rjet;
00165     double _JetPtCut;
00166     double _JetEtaCut;
00167     double _Lep1PtCut;
00168     double _Lep2PtCut;
00169     double _LepEtaCut;
00170     double _sumWeightSelected;
00171 
00172     //@{
00173     /// Histograms
00174     AIDA::IHistogram1D* _dStot;
00175     AIDA::IHistogram1D* _dSdET;
00176     AIDA::IHistogram1D* _dSdETA;
00177     AIDA::IHistogram1D* _dSdNJet;
00178     AIDA::IHistogram1D* _dSdNbJet;
00179     AIDA::IHistogram1D* _dSdZpT;
00180 
00181     //@}
00182 
00183   };
00184 
00185 
00186   // This global object acts as a hook for the plugin system
00187   AnalysisBuilder<CDF_2008_S8095620> plugin_CDF_2008_S8095620;
00188 
00189 }