CMS 3D CMS Logo

SoftLepton Class Reference

Description: CMSSW EDProducer for soft lepton b tagging. More...

#include <RecoBTag/SoftLepton/plugin/SoftLepton.h>

Inheritance diagram for SoftLepton:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 SoftLepton (const edm::ParameterSet &iConfig)
reco::SoftLeptonTagInfo tag (const edm::RefToBase< reco::Jet > &jet, const reco::TrackRefVector &tracks, const std::vector< edm::RefToBase< reco::Track > > &leptons, const reco::Vertex &primaryVertex) const
 ~SoftLepton ()

Protected Member Functions

GlobalVector refineJetAxis (const edm::RefToBase< reco::Jet > &jet, const reco::TrackRefVector &tracks, const edm::RefToBase< reco::Track > &exclude=edm::RefToBase< reco::Track >()) const

Static Protected Member Functions

static double boostedPPar (const math::XYZVector &vector, const math::XYZVector &axis)
static double relativeEta (const math::XYZVector &vector, const math::XYZVector &axis)

Private Member Functions

virtual void beginJob (const edm::EventSetup &setup)
virtual void endJob (void)
virtual void produce (edm::Event &event, const edm::EventSetup &setup)

Private Attributes

double m_chi2Cut
double m_deltaRCut
const edm::InputTag m_jets
const edm::InputTag m_leptons
const edm::InputTag m_primaryVertex
double m_qualityCut
unsigned int m_refineJetAxis
const TransientTrackBuilderm_transientTrackBuilder

Static Private Attributes

static const reco::Vertex s_nominalBeamSpot


Detailed Description

Description: CMSSW EDProducer for soft lepton b tagging.

Description: CMSSW EDProducer wrapper for sot lepton b tagging.

Implementation:

Implementation: The actual tagging is performed by SoftLeptonAlgorithm.

Definition at line 43 of file SoftLepton.h.


Constructor & Destructor Documentation

SoftLepton::SoftLepton ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 86 of file SoftLepton.cc.

00086                                                       :
00087   m_jets(          iConfig.getParameter<edm::InputTag>( "jets" ) ),
00088   m_primaryVertex( iConfig.getParameter<edm::InputTag>( "primaryVertex" ) ),
00089   m_leptons(       iConfig.getParameter<edm::InputTag>( "leptons" ) ),
00090   m_transientTrackBuilder( NULL ),
00091   m_refineJetAxis( iConfig.getParameter<unsigned int>( "refineJetAxis" ) ),
00092   m_deltaRCut(     iConfig.getParameter<double>( "leptonDeltaRCut" ) ),
00093   m_chi2Cut(       iConfig.getParameter<double>( "leptonChi2Cut" ) ),
00094   m_qualityCut(    iConfig.getParameter<double>( "leptonQualityCut" ) )
00095 {
00096   produces<reco::SoftLeptonTagInfoCollection>();
00097 }

SoftLepton::~SoftLepton ( void   ) 

Definition at line 100 of file SoftLepton.cc.

00100                             {
00101 }


Member Function Documentation

void SoftLepton::beginJob ( const edm::EventSetup setup  )  [private, virtual]

Reimplemented from edm::EDProducer.

Definition at line 228 of file SoftLepton.cc.

00228                                                 {
00229 }

double SoftLepton::boostedPPar ( const math::XYZVector vector,
const math::XYZVector axis 
) [static, protected]

Definition at line 357 of file SoftLepton.cc.

References metsig::jet.

00357                                                                                      {
00358   static const double lepton_mass = 0.00;       // assume a massless (ultrarelativistic) lepton
00359   static const double jet_mass    = 5.279;      // use B±/B0 mass as the jet rest mass [PDG 2007 updates]
00360   ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > lepton(vector.Dot(axis) / axis.r(), Perp(vector, axis), 0., lepton_mass);
00361   ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > jet( axis.r(), 0., 0., jet_mass );
00362   ROOT::Math::BoostX boost( -jet.Beta() );
00363   return boost(lepton).x();
00364 } 

void SoftLepton::endJob ( void   )  [private, virtual]

Reimplemented from edm::EDProducer.

Definition at line 233 of file SoftLepton.cc.

00233                        {
00234 }

void SoftLepton::produce ( edm::Event event,
const edm::EventSetup setup 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 105 of file SoftLepton.cc.

References and, metsig::electron, edm::EventSetup::get(), i, edm::Handle< T >::isValid(), pfTauBenchmarkGeneric_cfi::jets, edm::InputTag::label(), EgammaValidation_Zee_cff::leptons, m_jets, m_leptons, m_primaryVertex, m_qualityCut, m_transientTrackBuilder, metsig::muon, edm::errors::NotFound, edm::ESHandle< T >::product(), HLT_VtxMuL3::result, s_nominalBeamSpot, size, tag(), and tracks.

00105                                                                  {
00106 
00107   // grab a TransientTrack helper from the Event Setup
00108   edm::ESHandle<TransientTrackBuilder> builder;
00109   setup.get<TransientTrackRecord>().get( "TransientTrackBuilder", builder );
00110   m_transientTrackBuilder = builder.product();
00111 
00112   // input objects
00113 
00114   // input jets (and possibly tracks)
00115   ProductID jets_id;
00116   std::vector<edm::RefToBase<reco::Jet> > jets;
00117   std::vector<reco::TrackRefVector>       tracks;
00118   do { {
00119     // look for a JetTracksAssociationCollection
00120     edm::Handle<reco::JetTracksAssociationCollection> h_jtas;
00121     event.getByLabel(m_jets, h_jtas);
00122     if (h_jtas.isValid()) {
00123       unsigned int size = h_jtas->size();
00124       jets.resize(size);
00125       tracks.resize(size);
00126       for (unsigned int i = 0; i < size; ++i) {
00127         jets[i]   = (*h_jtas)[i].first;
00128         tracks[i] = (*h_jtas)[i].second;
00129       }
00130       break;
00131     }
00132   } { // else...
00133     // look for a View<Jet>
00134     edm::Handle<edm::View<reco::Jet> > h_jets;
00135     event.getByLabel(m_jets, h_jets);
00136     if (h_jets.isValid()) {
00137       unsigned int size = h_jets->size();
00138       jets.resize(size);
00139       tracks.resize(size);
00140       for (unsigned int i = 0; i < h_jets->size(); i++)
00141         jets[i] = h_jets->refAt(i);
00142       break;
00143     }
00144   } { // else...
00145     throw edm::Exception(edm::errors::NotFound) << "Object " << m_jets << " of type among (\"reco::JetTracksAssociationCollection\", \"edm::View<reco::Jet>\") not found";
00146   } } while (false);
00147   
00148   // input primary vetex (optional, can be "nominal" or "beamspot")
00149   reco::Vertex vertex;
00150   Handle<reco::VertexCollection> h_primaryVertex;
00151   if (m_primaryVertex.label() == "nominal") {
00152     vertex = s_nominalBeamSpot;
00153   } else
00154   if (m_primaryVertex.label() == "beamspot") {
00155     edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00156     event.getByType(recoBeamSpotHandle);
00157     vertex = reco::Vertex(recoBeamSpotHandle->position(), recoBeamSpotHandle->covariance3D(), 1, 1, 0);
00158   } else {
00159     event.getByLabel(m_primaryVertex, h_primaryVertex);
00160     if (h_primaryVertex->size()) {
00161       vertex = h_primaryVertex->front();
00162     } else {
00163       // fall back to nominal beam spot
00164       vertex = s_nominalBeamSpot;
00165     }
00166   }
00167 
00168   // input leptons (can be of different types)
00169   ProductID leptons_id;
00170   std::vector<edm::RefToBase<reco::Track> > leptons;
00171   // try to access the input collection as a collection of Electrons, Muons or Tracks
00172   do { {
00173     // look for vector<Electron>
00174     Handle<reco::ElectronCollection> h_electrons;
00175     event.getByLabel(m_leptons, h_electrons);
00176     if (h_electrons.isValid()) {
00177       for (reco::ElectronCollection::const_iterator electron = h_electrons->begin(); electron != h_electrons->end(); ++electron)
00178         leptons.push_back(edm::RefToBase<reco::Track>( electron->track() ));
00179       break;
00180     }
00181   } { // else
00182     // look for vector<GsfElectron>
00183     Handle<reco::GsfElectronCollection> h_electrons;
00184     event.getByLabel(m_leptons, h_electrons);
00185     if (h_electrons.isValid()) {
00186       for (reco::GsfElectronCollection::const_iterator electron = h_electrons->begin(); electron != h_electrons->end(); ++electron)
00187         leptons.push_back(edm::RefToBase<reco::Track>( electron->gsfTrack() ));
00188       break;
00189     }
00190   } { // else
00191     // look for vetor<Muon>
00192     Handle<reco::MuonCollection> h_muons;
00193     event.getByLabel(m_leptons, h_muons);
00194     if (h_muons.isValid()) {
00195       for (reco::MuonCollection::const_iterator muon = h_muons->begin(); muon != h_muons->end(); ++muon) {
00196         if (muon->isGlobalMuon())
00197           leptons.push_back(edm::RefToBase<reco::Track>( muon->globalTrack() ));
00198         else 
00199         if (muon->isTrackerMuon() and muon->caloCompatibility() > m_qualityCut)
00200           leptons.push_back(edm::RefToBase<reco::Track>( muon->innerTrack() ));
00201       }
00202       break;
00203     }
00204   } { // else
00205     // look for edm::View<Track> 
00206     Handle<edm::View<reco::Track> > h_tracks;
00207     event.getByLabel(m_leptons, h_tracks);
00208     if (h_tracks.isValid()) {
00209       for (unsigned int i = 0; i < h_tracks->size(); i++)
00210         leptons.push_back(h_tracks->refAt(i));
00211       break;
00212     }
00213   } { // else
00214     throw edm::Exception(edm::errors::NotFound) << "Object " << m_leptons << " of type among (\"reco::ElectronCollection\", \"reco::GsfElectronCollection\", \"reco::MuonCollection\", \"edm::View<reco::Track>\") not found";
00215   } } while (false);
00216 
00217   // output collections
00218   std::auto_ptr<reco::SoftLeptonTagInfoCollection> outputCollection(  new reco::SoftLeptonTagInfoCollection() );
00219   for (unsigned int i = 0; i < jets.size(); ++i) {
00220     reco::SoftLeptonTagInfo result = tag( jets[i], tracks[i], leptons, vertex );
00221     outputCollection->push_back( result );
00222   }
00223   event.put( outputCollection );
00224 }

GlobalVector SoftLepton::refineJetAxis ( const edm::RefToBase< reco::Jet > &  jet,
const reco::TrackRefVector tracks,
const edm::RefToBase< reco::Track > &  exclude = edm::RefToBase<reco::Track>() 
) const [protected]

Definition at line 276 of file SoftLepton.cc.

References and, AXIS_CHARGED_AVERAGE, AXIS_CHARGED_AVERAGE_NOLEPTON, AXIS_CHARGED_SUM, AXIS_CHARGED_SUM_NOLEPTON, edm::RefVector< C, T, F >::begin(), edm::RefVector< C, T, F >::end(), reco::TrackBase::eta(), edm::RefToBase< T >::isNonnull(), m_refineJetAxis, reco::TrackBase::momentum(), muonGeometry::perp(), reco::TrackBase::phi(), reco::TrackBase::pt(), sum(), and track.

00280         {
00281   math::XYZVector axis = jet->momentum();
00282 
00283   if (m_refineJetAxis == AXIS_CHARGED_AVERAGE or
00284       m_refineJetAxis == AXIS_CHARGED_AVERAGE_NOLEPTON) {
00285 
00286     double sum_pT        = 0.;
00287     double sum_eta_by_pT = 0.;
00288     double sum_phi_by_pT = 0.;
00289 
00290     double perp;
00291     double phi_rel;
00292     double eta_rel;
00293 
00294     // refine jet eta and phi with charged tracks measurements, if available
00295     for (reco::TrackRefVector::const_iterator track_it = tracks.begin(); track_it != tracks.end(); ++track_it ) {
00296       const reco::Track & track = **track_it;
00297 
00298       perp = track.pt();
00299       eta_rel = (double) track.eta() - axis.eta();
00300       phi_rel = (double) track.phi() - axis.phi();
00301       while (phi_rel < -M_PI) phi_rel += 2*M_PI;
00302       while (phi_rel >  M_PI) phi_rel -= 2*M_PI;
00303 
00304       sum_pT        += perp;
00305       sum_phi_by_pT += perp * phi_rel;
00306       sum_eta_by_pT += perp * eta_rel;
00307     }
00308 
00309     // "remove" excluded track
00310     if (m_refineJetAxis == AXIS_CHARGED_AVERAGE_NOLEPTON and exclude.isNonnull()) {
00311       const reco::Track & track = *exclude;
00312 
00313       perp = track.pt();
00314       eta_rel = (double) track.eta() - axis.eta();
00315       phi_rel = (double) track.phi() - axis.phi();
00316       while (phi_rel < -M_PI) phi_rel += 2*M_PI;
00317       while (phi_rel >  M_PI) phi_rel -= 2*M_PI;
00318 
00319       sum_pT        -= perp;
00320       sum_phi_by_pT -= perp * phi_rel;
00321       sum_eta_by_pT -= perp * eta_rel;
00322     }
00323 
00324     if (sum_pT > 1.)    // avoid the case of only the lepton-track with small rounding errors
00325       axis = math::RhoEtaPhiVector( axis.rho(), axis.eta() + sum_eta_by_pT / sum_pT, axis.phi() + sum_phi_by_pT / sum_pT);
00326     
00327   } else if (m_refineJetAxis == AXIS_CHARGED_SUM or
00328              m_refineJetAxis == AXIS_CHARGED_SUM_NOLEPTON) {
00329     math::XYZVector sum;
00330 
00331     // recalculate the jet direction as the sum of charget tracks momenta
00332     for (reco::TrackRefVector::const_iterator track_it = tracks.begin(); track_it != tracks.end(); ++track_it ) {
00333       const reco::Track & track = **track_it;
00334       sum += track.momentum();
00335     }
00336 
00337     // "remove" excluded track
00338     if (m_refineJetAxis == AXIS_CHARGED_SUM_NOLEPTON and exclude.isNonnull()) {
00339       const reco::Track & track = *exclude;
00340       sum -= track.momentum();
00341     }
00342 
00343     if (sum.R() > 1.) // avoid the case of only the lepton-track with small rounding errors
00344       axis = sum;
00345   }
00346   
00347   return GlobalVector(axis.x(), axis.y(), axis.z());
00348 }

double SoftLepton::relativeEta ( const math::XYZVector vector,
const math::XYZVector axis 
) [static, protected]

Definition at line 350 of file SoftLepton.cc.

References dot(), funct::log(), and muonGeometry::mag().

00350                                                                                      {
00351   double mag = vector.r() * axis.r();
00352   double dot = vector.Dot(axis); 
00353   return -log((mag - dot)/(mag + dot)) / 2;
00354 }

reco::SoftLeptonTagInfo SoftLepton::tag ( const edm::RefToBase< reco::Jet > &  jet,
const reco::TrackRefVector tracks,
const std::vector< edm::RefToBase< reco::Track > > &  leptons,
const reco::Vertex primaryVertex 
) const

Referenced by produce().


Member Data Documentation

double SoftLepton::m_chi2Cut [private]

Definition at line 89 of file SoftLepton.h.

double SoftLepton::m_deltaRCut [private]

Definition at line 88 of file SoftLepton.h.

const edm::InputTag SoftLepton::m_jets [private]

Definition at line 79 of file SoftLepton.h.

Referenced by produce().

const edm::InputTag SoftLepton::m_leptons [private]

Definition at line 81 of file SoftLepton.h.

Referenced by produce().

const edm::InputTag SoftLepton::m_primaryVertex [private]

Definition at line 80 of file SoftLepton.h.

Referenced by produce().

double SoftLepton::m_qualityCut [private]

Definition at line 90 of file SoftLepton.h.

Referenced by produce().

unsigned int SoftLepton::m_refineJetAxis [private]

Definition at line 87 of file SoftLepton.h.

Referenced by refineJetAxis().

const TransientTrackBuilder* SoftLepton::m_transientTrackBuilder [private]

Definition at line 84 of file SoftLepton.h.

Referenced by produce().

const reco::Vertex SoftLepton::s_nominalBeamSpot [static, private]

Definition at line 93 of file SoftLepton.h.

Referenced by produce().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:32:40 2009 for CMSSW by  doxygen 1.5.4