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/MergedFinalState.hh" 00006 #include "Rivet/Projections/DressedLeptons.hh" 00007 #include "Rivet/Projections/VetoedFinalState.hh" 00008 #include "Rivet/Projections/Beam.hh" 00009 00010 namespace Rivet { 00011 00012 00013 WFinder::WFinder(const FinalState& inputfs, 00014 const Cut& leptoncuts, 00015 PdgId pid, 00016 double minmass, double maxmass, 00017 double missingET, 00018 double dRmax, 00019 ClusterPhotons clusterPhotons, 00020 PhotonTracking trackPhotons, 00021 MassWindow masstype, 00022 double masstarget) { 00023 setName("WFinder"); 00024 00025 _minmass = minmass; 00026 _maxmass = maxmass; 00027 _masstarget = masstarget; 00028 _pid = pid; 00029 _trackPhotons = trackPhotons; 00030 _useTransverseMass = (masstype == TRANSMASS); 00031 00032 // Check that the arguments are legal 00033 assert(abs(_pid) == PID::ELECTRON || abs(_pid) == PID::MUON); 00034 _nu_pid = abs(_pid) + 1; 00035 assert(abs(_nu_pid) == PID::NU_E || abs(_nu_pid) == PID::NU_MU); 00036 00037 // Lepton clusters 00038 IdentifiedFinalState bareleptons(inputfs); 00039 bareleptons.acceptIdPair(pid); 00040 const bool doClustering = (clusterPhotons != NOCLUSTER); 00041 const bool useDecayPhotons = (clusterPhotons == CLUSTERALL); 00042 DressedLeptons leptons(inputfs, bareleptons, dRmax, leptoncuts, doClustering, useDecayPhotons); 00043 addProjection(leptons, "DressedLeptons"); 00044 00045 // Add MissingMomentum proj to calc MET 00046 MissingMomentum vismom(inputfs); 00047 addProjection(vismom, "MissingET"); 00048 // Set ETmiss cut 00049 _etMiss = missingET; 00050 00051 VetoedFinalState remainingFS; 00052 remainingFS.addVetoOnThisFinalState(*this); 00053 addProjection(remainingFS, "RFS"); 00054 } 00055 00056 00057 ///////////////////////////////////////////////////// 00058 00059 00060 const VetoedFinalState& WFinder::remainingFinalState() const { 00061 return getProjection<VetoedFinalState>("RFS"); 00062 } 00063 00064 00065 const MissingMomentum& WFinder::missingMom() const { 00066 return getProjection<MissingMomentum>("MissingET"); 00067 } 00068 00069 00070 int WFinder::compare(const Projection& p) const { 00071 PCmp dlcmp = mkNamedPCmp(p, "DressedLeptons"); 00072 if (dlcmp != EQUIVALENT) return dlcmp; 00073 00074 const WFinder& other = dynamic_cast<const WFinder&>(p); 00075 return (cmp(_minmass, other._minmass) || cmp(_maxmass, other._maxmass) || 00076 cmp(_useTransverseMass, other._useTransverseMass) || 00077 cmp(_etMiss, other._etMiss) || 00078 cmp(_pid, other._pid) || cmp(_trackPhotons, other._trackPhotons)); 00079 } 00080 00081 00082 void WFinder::project(const Event& e) { 00083 clear(); 00084 00085 // Beam beam; 00086 // beam.project(e); 00087 // const double sqrtS = beam.sqrtS(); 00088 00089 // Check missing ET 00090 const MissingMomentum& missmom = applyProjection<MissingMomentum>(e, "MissingET"); 00091 //const FourMomentum& pmiss = FourMomentum(sqrtS,0,0,0) + missmom.missingMomentum(); 00092 00093 const double met = missmom.vectorEt().mod(); 00094 MSG_TRACE("MET = " << met/GeV << " GeV vs. required = " << _etMiss/GeV << " GeV"); 00095 if (met < _etMiss) { 00096 MSG_DEBUG("Not enough missing ET: " << met/GeV << " GeV vs. " << _etMiss/GeV << " GeV"); 00097 return; 00098 } 00099 00100 const DressedLeptons& leptons = applyProjection<DressedLeptons>(e, "DressedLeptons"); 00101 if ( leptons.dressedLeptons().empty() ) { 00102 MSG_DEBUG("No dressed leptons."); 00103 return; 00104 } 00105 00106 MSG_DEBUG("Found at least one dressed lepton: " << leptons.dressedLeptons()[0].momentum() ); 00107 const FourMomentum pmiss = missmom.missingMomentum(0*GeV); 00108 MSG_DEBUG("Found missing 4-momentum: " << pmiss); 00109 00110 // Make and register an invariant mass final state for the W decay leptons 00111 vector<pair<PdgId, PdgId> > l_nu_ids; 00112 l_nu_ids += make_pair(abs(_pid), -_nu_pid); 00113 l_nu_ids += make_pair(-abs(_pid), _nu_pid); 00114 InvMassFinalState imfs(l_nu_ids, _minmass, _maxmass, _masstarget); 00115 imfs.useTransverseMass(_useTransverseMass); 00116 Particles tmp; 00117 tmp.insert(tmp.end(), leptons.dressedLeptons().begin(), leptons.dressedLeptons().end()); 00118 tmp.push_back(Particle( _nu_pid, pmiss)); // fake neutrino from Et miss vector 00119 tmp.push_back(Particle(-_nu_pid, pmiss)); // fake antineutrino from Et miss vector 00120 imfs.calc(tmp); 00121 00122 if (imfs.particlePairs().size() < 1) return; 00123 00124 ParticlePair Wconstituents(imfs.particlePairs()[0]); 00125 Particle p1(Wconstituents.first), p2(Wconstituents.second); 00126 00127 if (threeCharge(p1) == 0) { 00128 _constituentLeptons += p2; 00129 _constituentNeutrinos += p1; 00130 } else { 00131 _constituentLeptons += p1; 00132 _constituentNeutrinos += p2; 00133 } 00134 00135 FourMomentum pW = p1.momentum() + p2.momentum(); 00136 const int w3charge = threeCharge(p1) + threeCharge(p2); 00137 assert(abs(w3charge) == 3); 00138 const int wcharge = w3charge/3; 00139 00140 stringstream msg; 00141 string wsign = (wcharge == 1) ? "+" : "-"; 00142 string wstr = "W" + wsign; 00143 msg << wstr << " " << pW << " reconstructed from: " << "\n" 00144 << " " << p1.momentum() << " " << p1.pid() << "\n" 00145 << " + " << p2.momentum() << " " << p2.pid(); 00146 MSG_DEBUG(msg.str()); 00147 00148 // Make W Particle and insert into particles list 00149 const PdgId wpid = (wcharge == 1) ? PID::WPLUSBOSON : PID::WMINUSBOSON; 00150 _bosons.push_back(Particle(wpid, pW)); 00151 00152 // Find the DressedLeptons which survived the IMFS cut such that we can 00153 // extract their original particles 00154 00155 // TODO: do we need to add all used invisibles to _theParticles ? 00156 00157 foreach (const Particle& p, _constituentLeptons) { 00158 foreach (const DressedLepton& l, leptons.dressedLeptons()) { 00159 if (p.pid() == l.pid() && p.momentum() == l.momentum()) { 00160 _theParticles.push_back(l.constituentLepton()); 00161 if (_trackPhotons) { 00162 _theParticles.insert(_theParticles.end(), l.constituentPhotons().begin(), l.constituentPhotons().end()); 00163 } 00164 } 00165 } 00166 } 00167 } 00168 00169 00170 } Generated on Wed Oct 7 2015 12:09:15 for The Rivet MC analysis system by ![]() |