rivet is hosted by Hepforge, IPPP Durham
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 
00013   WFinder::WFinder(const FinalState& inputfs,
00014                    const Cut& fsCut,
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     // Don't make pT or eta cuts on the neutrino
00038     IdentifiedFinalState neutrinos(inputfs);
00039     neutrinos.acceptNeutrinos();
00040     addProjection(neutrinos, "Neutrinos");
00041 
00042     // Lepton clusters
00043     IdentifiedFinalState bareleptons(inputfs);
00044     bareleptons.acceptIdPair(pid);
00045     const bool doClustering = (clusterPhotons != NOCLUSTER);
00046     const bool useDecayPhotons = (clusterPhotons == CLUSTERALL);
00047     DressedLeptons leptons(inputfs, bareleptons, dRmax, fsCut, doClustering, useDecayPhotons);
00048     addProjection(leptons, "DressedLeptons");
00049 
00050     // Add MissingMomentum proj to calc MET
00051     MissingMomentum vismom(inputfs);
00052     addProjection(vismom, "MissingET");
00053     // Set ETmiss
00054     _etMiss = missingET;
00055 
00056     VetoedFinalState remainingFS;
00057     remainingFS.addVetoOnThisFinalState(*this);
00058     addProjection(remainingFS, "RFS");
00059   }
00060 
00061 
00062   /////////////////////////////////////////////////////
00063 
00064 
00065   const FinalState& WFinder::remainingFinalState() const {
00066     return getProjection<FinalState>("RFS");
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     const DressedLeptons& leptons = applyProjection<DressedLeptons>(e, "DressedLeptons");
00086     const FinalState& neutrinos = applyProjection<FinalState>(e, "Neutrinos");
00087 
00088     // Make and register an invariant mass final state for the W decay leptons
00089     vector<pair<PdgId, PdgId> > l_nu_ids;
00090     l_nu_ids += make_pair(abs(_pid), -abs(_nu_pid));
00091     l_nu_ids += make_pair(-abs(_pid), abs(_nu_pid));
00092     InvMassFinalState imfs(l_nu_ids, _minmass, _maxmass, _masstarget);
00093     imfs.useTransverseMass(_useTransverseMass);
00094     Particles tmp;
00095     tmp.insert(tmp.end(), leptons.dressedLeptons().begin(), leptons.dressedLeptons().end());
00096     tmp.insert(tmp.end(), neutrinos.particles().begin(), neutrinos.particles().end());
00097     imfs.calc(tmp);
00098 
00099     if (imfs.particlePairs().size() < 1) return;
00100 
00101     ParticlePair Wconstituents(imfs.particlePairs()[0]);
00102     Particle p1(Wconstituents.first), p2(Wconstituents.second);
00103 
00104     if (threeCharge(p1) == 0) {
00105       _constituentLeptons += p2;
00106       _constituentNeutrinos += p1;
00107     } else {
00108       _constituentLeptons += p1;
00109       _constituentNeutrinos += p2;
00110     }
00111 
00112     FourMomentum pW = p1.momentum() + p2.momentum();
00113     const int w3charge = threeCharge(p1) + threeCharge(p2);
00114     assert(abs(w3charge) == 3);
00115     const int wcharge = w3charge/3;
00116 
00117     stringstream msg;
00118     string wsign = (wcharge == 1) ? "+" : "-";
00119     string wstr = "W" + wsign;
00120     msg << wstr << " reconstructed from: " << "\n"
00121         << "   " << p1.momentum() << " " << p1.pid() << "\n"
00122         << " + " << p2.momentum() << " " << p2.pid();
00123     MSG_DEBUG(msg.str());
00124 
00125     // Check missing ET
00126     const MissingMomentum& vismom = applyProjection<MissingMomentum>(e, "MissingET");
00127     const double met = vismom.vectorEt().mod();
00128     MSG_TRACE("MET = " << met/GeV << " GeV vs. required = " << _etMiss/GeV << " GeV");
00129     if (met < _etMiss) {
00130       MSG_DEBUG("Not enough missing ET: " << met/GeV << " GeV vs. " << _etMiss/GeV << " GeV");
00131       return;
00132     }
00133 
00134     // Make W Particle and insert into particles list
00135     const PdgId wpid = (wcharge == 1) ? PID::WPLUSBOSON : PID::WMINUSBOSON;
00136     _bosons.push_back(Particle(wpid, pW));
00137 
00138     // Find the DressedLeptons and neutrinos which survived the IMFS cut such that we can
00139     // extract their original particles
00140     foreach (const Particle& p, _constituentNeutrinos) {
00141       _theParticles.push_back(p);
00142     }
00143     foreach (const Particle& p, _constituentLeptons) {
00144       foreach (const DressedLepton& l, leptons.dressedLeptons()) {
00145         if (p.pid() == l.pid() && p.momentum() == l.momentum()) {
00146           _theParticles.push_back(l.constituentLepton());
00147           if (_trackPhotons) {
00148             _theParticles.insert(_theParticles.end(), l.constituentPhotons().begin(), l.constituentPhotons().end());
00149           }
00150         }
00151       }
00152     }
00153   }
00154 
00155 
00156 }