CDF_2006_S6653332.cc
Go to the documentation of this file.00001
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/RivetAIDA.hh"
00004 #include "Rivet/Tools/Logging.hh"
00005 #include "Rivet/Tools/ParticleIdUtils.hh"
00006 #include "Rivet/Projections/FinalState.hh"
00007 #include "Rivet/Projections/VetoedFinalState.hh"
00008 #include "Rivet/Projections/ChargedFinalState.hh"
00009 #include "Rivet/Projections/InvMassFinalState.hh"
00010 #include "Rivet/Projections/FastJets.hh"
00011 #include "Rivet/Projections/ChargedLeptons.hh"
00012
00013 namespace Rivet {
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 class CDF_2006_S6653332 : public Analysis {
00024 public:
00025
00026
00027 CDF_2006_S6653332()
00028 : Analysis("CDF_2006_S6653332"),
00029 _Rjet(0.7), _JetPtCut(20.), _JetEtaCut(1.5), _Lep1PtCut(18.), _Lep2PtCut(10.), _LepEtaCut(1.1),
00030 _sumWeightsWithZ(0.0), _sumWeightsWithZJet(0.0)
00031 {
00032 setBeams(PROTON, ANTIPROTON);
00033 setNeedsCrossSection(true);
00034 }
00035
00036
00037
00038
00039
00040 void init() {
00041 const FinalState fs(-3.6, 3.6);
00042 addProjection(fs, "FS");
00043
00044
00045
00046 vector<pair<PdgId,PdgId> > vids;
00047 vids.push_back(make_pair(ELECTRON, POSITRON));
00048 vids.push_back(make_pair(MUON, ANTIMUON));
00049 FinalState fs2(-3.6, 3.6);
00050 InvMassFinalState invfs(fs2, vids, 66*GeV, 116*GeV);
00051 addProjection(invfs, "INVFS");
00052
00053
00054 VetoedFinalState vfs(fs);
00055 vfs.addVetoOnThisFinalState(invfs);
00056 addProjection(vfs, "VFS");
00057 addProjection(FastJets(vfs, FastJets::CDFMIDPOINT, 0.7), "Jets");
00058
00059
00060 _sigmaBJet = bookHistogram1D(1, 1, 1);
00061 _ratioBJetToZ = bookHistogram1D(2, 1, 1);
00062 _ratioBJetToJet = bookHistogram1D(3, 1, 1);
00063 }
00064
00065
00066
00067 void analyze(const Event& event) {
00068
00069
00070 const InvMassFinalState& invMassFinalState = applyProjection<InvMassFinalState>(event, "INVFS");
00071 const ParticleVector& ZDecayProducts = invMassFinalState.particles();
00072
00073
00074 if (ZDecayProducts.size() < 2) vetoEvent;
00075
00076 double Lep1Pt = ZDecayProducts[0].momentum().perp();
00077 double Lep2Pt = ZDecayProducts[1].momentum().perp();
00078 double Lep1Eta = fabs(ZDecayProducts[0].momentum().rapidity());
00079 double Lep2Eta = fabs(ZDecayProducts[1].momentum().rapidity());
00080
00081 if (Lep1Eta > _LepEtaCut && Lep2Eta > _LepEtaCut) vetoEvent;
00082
00083 if (abs(ZDecayProducts[0].pdgId())==13 && (Lep1Eta > 1.0 && Lep2Eta > 1.)) {
00084 vetoEvent;
00085 }
00086 if (Lep1Pt < _Lep1PtCut && Lep2Pt < _Lep1PtCut) vetoEvent;
00087
00088
00089
00090 _sumWeightsWithZ += event.weight();
00091
00092 FourMomentum Zmom = ZDecayProducts[0].momentum() + ZDecayProducts[1].momentum();
00093
00094
00095
00096 ParticleVector bquarks;
00097
00098 for (GenEvent::particle_const_iterator p = event.genEvent().particles_begin();
00099 p != event.genEvent().particles_end(); ++p) {
00100 if ( fabs((*p)->pdg_id()) == BQUARK ) {
00101 bquarks.push_back(Particle(**p));
00102 }
00103 }
00104
00105
00106 const FastJets& jetpro = applyProjection<FastJets>(event, "Jets");
00107 getLog() << Log::DEBUG << "Jet multiplicity before any pT cut = " << jetpro.size() << endl;
00108
00109 const PseudoJets& jets = jetpro.pseudoJetsByPt();
00110 getLog() << Log::DEBUG << "jetlist size = " << jets.size() << endl;
00111
00112 int numBJet = 0;
00113 int numJet = 0;
00114
00115
00116 for (PseudoJets::const_iterator jt = jets.begin(); jt != jets.end(); ++jt) {
00117
00118 if (jt->perp() > _JetPtCut && fabs(jt->rapidity()) <= _JetEtaCut) {
00119 ++numJet;
00120
00121
00122
00123 bool bjet = false;
00124 foreach (const Particle& bquark, bquarks) {
00125 if (deltaR(jt->rapidity(), jt->phi(), bquark.momentum().rapidity(),bquark.momentum().azimuthalAngle()) <= _Rjet) {
00126 bjet = true;
00127 break;
00128 }
00129 }
00130 if (bjet) {
00131 numBJet++;
00132 }
00133 }
00134 }
00135
00136 if (numJet > 0) _sumWeightsWithZJet += event.weight();
00137 if (numBJet > 0) {
00138 _sigmaBJet->fill(1960.0,event.weight());
00139 _ratioBJetToZ->fill(1960.0,event.weight());
00140 _ratioBJetToJet->fill(1960.0,event.weight());
00141 }
00142
00143 }
00144
00145
00146
00147 void finalize() {
00148 getLog() << Log::DEBUG << "Total sum of weights = " << sumOfWeights() << endl;
00149 getLog() << Log::DEBUG << "Sum of weights for Z production in mass range = " << _sumWeightsWithZ << endl;
00150 getLog() << Log::DEBUG << "Sum of weights for Z+jet production in mass range = " << _sumWeightsWithZJet << endl;
00151
00152 _sigmaBJet->scale(crossSection()/sumOfWeights());
00153 _ratioBJetToZ->scale(1.0/_sumWeightsWithZ);
00154 _ratioBJetToJet->scale(1.0/_sumWeightsWithZJet);
00155 }
00156
00157
00158
00159
00160 private:
00161
00162
00163
00164
00165 double _Rjet;
00166 double _JetPtCut;
00167 double _JetEtaCut;
00168 double _Lep1PtCut;
00169 double _Lep2PtCut;
00170 double _LepEtaCut;
00171
00172 double _sumWeightsWithZ;
00173 double _sumWeightsWithZJet;
00174
00175
00176
00177
00178
00179 AIDA::IHistogram1D* _sigmaBJet;
00180 AIDA::IHistogram1D* _ratioBJetToZ;
00181 AIDA::IHistogram1D* _ratioBJetToJet;
00182
00183
00184 };
00185
00186
00187
00188 AnalysisBuilder<CDF_2006_S6653332> plugin_CDF_2006_S6653332;
00189
00190 }