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