WFinder.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Projections/WFinder.hh" 00003 #include "Rivet/Projections/InvMassFinalState.hh" 00004 #include "Rivet/Projections/MissingMomentum.hh" 00005 #include "Rivet/Projections/MergedFinalState.hh" 00006 #include "Rivet/Projections/DressedLeptons.hh" 00007 #include "Rivet/Projections/VetoedFinalState.hh" 00008 00009 namespace Rivet { 00010 00011 WFinder::WFinder(const FinalState& inputfs, 00012 Cut fsCut, 00013 PdgId pid, 00014 double minmass, double maxmass, 00015 double missingET, 00016 double dRmax, 00017 ClusterPhotons clusterPhotons, 00018 PhotonTracking trackPhotons, 00019 MassWindow masstype, 00020 double masstarget) { 00021 setName("WFinder"); 00022 00023 _minmass = minmass; 00024 _maxmass = maxmass; 00025 _masstarget = masstarget; 00026 _pid = pid; 00027 _trackPhotons = trackPhotons; 00028 _useTransverseMass = (masstype == MASS); 00029 00030 // Check that the arguments are legal 00031 assert(abs(_pid) == PID::ELECTRON || abs(_pid) == PID::MUON); 00032 _nu_pid = abs(_pid) + 1; 00033 assert(abs(_nu_pid) == PID::NU_E || abs(_nu_pid) == PID::NU_MU); 00034 00035 // Don't make pT or eta cuts on the neutrino 00036 IdentifiedFinalState neutrinos(inputfs); 00037 neutrinos.acceptNeutrinos(); 00038 addProjection(neutrinos, "Neutrinos"); 00039 00040 // Lepton clusters 00041 IdentifiedFinalState bareleptons(inputfs); 00042 bareleptons.acceptIdPair(pid); 00043 const bool doClustering = (clusterPhotons != NOCLUSTER); 00044 const bool useDecayPhotons = (clusterPhotons == CLUSTERALL); 00045 DressedLeptons leptons(inputfs, bareleptons, dRmax, doClustering, fsCut, useDecayPhotons); 00046 addProjection(leptons, "DressedLeptons"); 00047 00048 // Add MissingMomentum proj to calc MET 00049 MissingMomentum vismom(inputfs); 00050 addProjection(vismom, "MissingET"); 00051 // Set ETmiss 00052 _etMiss = missingET; 00053 00054 VetoedFinalState remainingFS; 00055 remainingFS.addVetoOnThisFinalState(*this); 00056 addProjection(remainingFS, "RFS"); 00057 } 00058 00059 00060 ///////////////////////////////////////////////////// 00061 00062 00063 const FinalState& WFinder::remainingFinalState() const { 00064 return getProjection<FinalState>("RFS"); 00065 } 00066 00067 00068 int WFinder::compare(const Projection& p) const { 00069 PCmp LCcmp = mkNamedPCmp(p, "DressedLeptons"); 00070 if (LCcmp != EQUIVALENT) return LCcmp; 00071 00072 const WFinder& other = dynamic_cast<const WFinder&>(p); 00073 return (cmp(_minmass, other._minmass) || cmp(_maxmass, other._maxmass) || 00074 cmp(_useTransverseMass, other._useTransverseMass) || 00075 cmp(_etMiss, other._etMiss) || 00076 cmp(_pid, other._pid) || cmp(_trackPhotons, other._trackPhotons)); 00077 } 00078 00079 00080 void WFinder::project(const Event& e) { 00081 clear(); 00082 00083 const DressedLeptons& leptons = applyProjection<DressedLeptons>(e, "DressedLeptons"); 00084 const FinalState& neutrinos = applyProjection<FinalState>(e, "Neutrinos"); 00085 00086 // Make and register an invariant mass final state for the W decay leptons 00087 vector<pair<PdgId, PdgId> > l_nu_ids; 00088 l_nu_ids += make_pair(abs(_pid), -abs(_nu_pid)); 00089 l_nu_ids += make_pair(-abs(_pid), abs(_nu_pid)); 00090 InvMassFinalState imfs(l_nu_ids, _minmass, _maxmass, _masstarget); 00091 imfs.useTransverseMass(_useTransverseMass); 00092 Particles tmp; 00093 tmp.insert(tmp.end(), leptons.clusteredLeptons().begin(), leptons.clusteredLeptons().end()); 00094 tmp.insert(tmp.end(), neutrinos.particles().begin(), neutrinos.particles().end()); 00095 imfs.calc(tmp); 00096 00097 if (imfs.particlePairs().size() < 1) return; 00098 00099 ParticlePair Wconstituents(imfs.particlePairs()[0]); 00100 Particle p1(Wconstituents.first), p2(Wconstituents.second); 00101 00102 if (PID::threeCharge(p1) == 0) { 00103 _constituentLeptons += p2; 00104 _constituentNeutrinos += p1; 00105 } else { 00106 _constituentLeptons += p1; 00107 _constituentNeutrinos += p2; 00108 } 00109 00110 FourMomentum pW = p1.momentum() + p2.momentum(); 00111 const int w3charge = PID::threeCharge(p1) + PID::threeCharge(p2); 00112 assert(abs(w3charge) == 3); 00113 const int wcharge = w3charge/3; 00114 00115 stringstream msg; 00116 string wsign = (wcharge == 1) ? "+" : "-"; 00117 string wstr = "W" + wsign; 00118 msg << wstr << " reconstructed from: " << "\n" 00119 << " " << p1.momentum() << " " << p1.pdgId() << "\n" 00120 << " + " << p2.momentum() << " " << p2.pdgId(); 00121 MSG_DEBUG(msg.str()); 00122 00123 // Check missing ET 00124 const MissingMomentum& vismom = applyProjection<MissingMomentum>(e, "MissingET"); 00125 /// @todo Restrict missing momentum eta range? Use vectorET()? 00126 if (vismom.scalarEt() < _etMiss) { 00127 MSG_DEBUG("Not enough missing ET: " << vismom.scalarEt()/GeV 00128 << " GeV vs. " << _etMiss/GeV << " GeV"); 00129 return; 00130 } 00131 00132 // Make W Particle and insert into particles list 00133 const PdgId wpid = (wcharge == 1) ? PID::WPLUSBOSON : PID::WMINUSBOSON; 00134 _bosons.push_back(Particle(wpid, pW)); 00135 00136 // Find the DressedLeptons and neutrinos which survived the IMFS cut such that we can 00137 // extract their original particles 00138 foreach (const Particle& p, _constituentNeutrinos) { 00139 _theParticles.push_back(p); 00140 } 00141 foreach (const Particle& p, _constituentLeptons) { 00142 foreach (const ClusteredLepton& l, leptons.clusteredLeptons()) { 00143 if (p.pdgId() == l.pdgId() && p.momentum() == l.momentum()) { 00144 _theParticles.push_back(l.constituentLepton()); 00145 if (_trackPhotons) { 00146 _theParticles.insert(_theParticles.end(), l.constituentPhotons().begin(), l.constituentPhotons().end()); 00147 } 00148 } 00149 } 00150 } 00151 } 00152 00153 00154 } Generated on Thu Feb 6 2014 17:38:47 for The Rivet MC analysis system by ![]() |