CMS 3D CMS Logo

Classes | Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Private Member Functions | Private Attributes | Static Private Attributes

SoftLepton Class Reference

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

Inheritance diagram for SoftLepton:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Classes

struct  TrackCompare

Public Types

typedef std::map< unsigned int,
float > 
LeptonIds
typedef std::map
< edm::RefToBase< reco::Track >
, LeptonIds, TrackCompare
Leptons

Public Member Functions

 SoftLepton (const edm::ParameterSet &iConfig)
reco::SoftLeptonTagInfo tag (const edm::RefToBase< reco::Jet > &jet, const reco::TrackRefVector &tracks, const Leptons &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 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_leptonCands
const edm::InputTag m_leptonId
const edm::InputTag m_leptons
muon::SelectionType m_muonSelection
const edm::InputTag m_primaryVertex
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.

Implementation:

Description: CMSSW EDProducer wrapper for sot lepton b tagging.

Implementation: The actual tagging is performed by SoftLeptonAlgorithm.

Definition at line 44 of file SoftLepton.h.


Member Typedef Documentation

typedef std::map<unsigned int, float> SoftLepton::LeptonIds

Definition at line 57 of file SoftLepton.h.

Definition at line 58 of file SoftLepton.h.


Constructor & Destructor Documentation

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

Definition at line 94 of file SoftLepton.cc.

                                                      :
  m_jets(          iConfig.getParameter<edm::InputTag>( "jets" ) ),
  m_primaryVertex( iConfig.getParameter<edm::InputTag>( "primaryVertex" ) ),
  m_leptons(       iConfig.getParameter<edm::InputTag>( "leptons" ) ),
  m_leptonCands(   iConfig.exists("leptonCands") ? iConfig.getParameter<edm::InputTag>( "leptonCands" ) : edm::InputTag() ),
  m_leptonId(      iConfig.exists("leptonId") ? iConfig.getParameter<edm::InputTag>( "leptonId" ) : edm::InputTag() ),
  m_transientTrackBuilder( 0 ),
  m_refineJetAxis( iConfig.getParameter<unsigned int>( "refineJetAxis" ) ),
  m_deltaRCut(     iConfig.getParameter<double>( "leptonDeltaRCut" ) ),
  m_chi2Cut(       iConfig.getParameter<double>( "leptonChi2Cut" ) ),
  m_muonSelection( (muon::SelectionType) iConfig.getParameter<unsigned int>( "muonSelection" ) )
{
  produces<reco::SoftLeptonTagInfoCollection>();
}
SoftLepton::~SoftLepton ( void  )

Definition at line 110 of file SoftLepton.cc.

                            {
}

Member Function Documentation

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

Definition at line 412 of file SoftLepton.cc.

References metsig::jet.

Referenced by tag().

                                                                                     {
  static const double lepton_mass = 0.00;       // assume a massless (ultrarelativistic) lepton
  static const double jet_mass    = 5.279;      // use B±/B0 mass as the jet rest mass [PDG 2007 updates]
  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > lepton(vector.Dot(axis) / axis.r(), Perp(vector, axis), 0., lepton_mass);
  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > jet( axis.r(), 0., 0., jet_mass );
  ROOT::Math::BoostX boost( -jet.Beta() );
  return boost(lepton).x();
} 
void SoftLepton::produce ( edm::Event event,
const edm::EventSetup setup 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 115 of file SoftLepton.cc.

References metsig::electron, Exception, edm::EventSetup::get(), i, errorMatrix2Lands_multiChannel::id, muon::isGoodMuon(), edm::HandleBase::isValid(), analyzePatCleaning_cfg::jets, EgammaValidation_Wenu_cff::leptons, m_jets, m_leptonCands, m_leptonId, m_leptons, m_muonSelection, m_primaryVertex, m_transientTrackBuilder, metsig::muon, edm::errors::NotFound, edm::ESHandle< T >::product(), query::result, s_nominalBeamSpot, findQualityFiles::size, tag(), and testEve_cfg::tracks.

                                                                 {

  // grab a TransientTrack helper from the Event Setup
  edm::ESHandle<TransientTrackBuilder> builder;
  setup.get<TransientTrackRecord>().get( "TransientTrackBuilder", builder );
  m_transientTrackBuilder = builder.product();

  // input objects

  // input jets (and possibly tracks)
  ProductID jets_id;
  std::vector<edm::RefToBase<reco::Jet> > jets;
  std::vector<reco::TrackRefVector>       tracks;
  do { {
    // look for a JetTracksAssociationCollection
    edm::Handle<reco::JetTracksAssociationCollection> h_jtas;
    event.getByLabel(m_jets, h_jtas);
    if (h_jtas.isValid()) {
      unsigned int size = h_jtas->size();
      jets.resize(size);
      tracks.resize(size);
      for (unsigned int i = 0; i < size; ++i) {
        jets[i]   = (*h_jtas)[i].first;
        tracks[i] = (*h_jtas)[i].second;
      }
      break;
    }
  } { // else...
    // look for a View<Jet>
    edm::Handle<edm::View<reco::Jet> > h_jets;
    event.getByLabel(m_jets, h_jets);
    if (h_jets.isValid()) {
      unsigned int size = h_jets->size();
      jets.resize(size);
      tracks.resize(size);
      for (unsigned int i = 0; i < h_jets->size(); i++)
        jets[i] = h_jets->refAt(i);
      break;
    }
  } { // else...
    throw edm::Exception(edm::errors::NotFound) << "Object " << m_jets << " of type among (\"reco::JetTracksAssociationCollection\", \"edm::View<reco::Jet>\") not found";
  } } while (false);
  
  // input primary vetex
  reco::Vertex vertex;
  Handle<reco::VertexCollection> h_primaryVertex;
  event.getByLabel(m_primaryVertex, h_primaryVertex);
  if (h_primaryVertex.isValid() and not h_primaryVertex->empty())
    vertex = h_primaryVertex->front();
  else
    // fall back to nominal beam spot
    vertex = s_nominalBeamSpot;

  // input leptons (can be of different types)
  Leptons leptons;

  Handle< ValueMap<float> > h_leptonCands;
  bool haveLeptonCands = !(m_leptonCands == edm::InputTag());
  if (haveLeptonCands)
    event.getByLabel(m_leptonCands, h_leptonCands);

  // try to access the input collection as a collection of GsfElectrons, Muons or Tracks

  unsigned int leptonId = SoftLeptonProperties::quality::leptonId;
  do { {
    // look for View<GsfElectron>
    Handle<GsfElectronView> h_electrons;
    event.getByLabel(m_leptons, h_electrons);

    if (h_electrons.isValid()) {
      leptonId = SoftLeptonProperties::quality::egammaElectronId;
      for (GsfElectronView::const_iterator electron = h_electrons->begin(); electron != h_electrons->end(); ++electron) {
        LeptonIds &id = leptons[reco::TrackBaseRef(electron->gsfTrack())];
        id[SoftLeptonProperties::quality::pfElectronId] = electron->mva();
        if (haveLeptonCands)
          id[SoftLeptonProperties::quality::btagElectronCands] = (*h_leptonCands)[h_electrons->refAt(electron - h_electrons->begin())];
      }
      break;
    }
  } { // else
    // look for View<Electron>
    // FIXME: is this obsolete?
    Handle<ElectronView> h_electrons;
    event.getByLabel(m_leptons, h_electrons);
    if (h_electrons.isValid()) {
      leptonId = SoftLeptonProperties::quality::egammaElectronId;
      for (ElectronView::const_iterator electron = h_electrons->begin(); electron != h_electrons->end(); ++electron) {
        LeptonIds &id = leptons[reco::TrackBaseRef(electron->track())];
        if (haveLeptonCands)
          id[SoftLeptonProperties::quality::btagElectronCands] = (*h_leptonCands)[h_electrons->refAt(electron - h_electrons->begin())];
      }
      break;
    }
  } { // else
    // look for PFElectrons
    // FIXME: is this obsolete?
    Handle<reco::PFCandidateCollection> h_electrons;
    event.getByLabel(m_leptons, h_electrons);
    if (h_electrons.isValid()) {
      leptonId = SoftLeptonProperties::quality::egammaElectronId;
      for (reco::PFCandidateCollection::const_iterator electron = h_electrons->begin(); electron != h_electrons->end(); ++electron) {
        LeptonIds *id;
        if (electron->gsfTrackRef().isNonnull())
          id = &leptons[reco::TrackBaseRef(electron->gsfTrackRef())];
        else if (electron->trackRef().isNonnull())
          id = &leptons[reco::TrackBaseRef(electron->trackRef())];
        else
          continue;
        (*id)[SoftLeptonProperties::quality::pfElectronId] = electron->mva_e_pi();
        if (haveLeptonCands)
          (*id)[SoftLeptonProperties::quality::btagElectronCands] = (*h_leptonCands)[reco::PFCandidateRef(h_electrons, electron - h_electrons->begin())];
      }
      break;
    }
  } { // else
    // look for View<Muon>
    Handle<MuonView> h_muons;
    event.getByLabel(m_leptons, h_muons);
    if (h_muons.isValid()) {
      for (MuonView::const_iterator muon = h_muons->begin(); muon != h_muons->end(); ++muon) {
        // FIXME -> turn this selection into a muonCands input?
        if (muon::isGoodMuon( *muon, m_muonSelection )) {
          LeptonIds *id;
          if (muon->globalTrack().isNonnull())
            id = &leptons[reco::TrackBaseRef(muon->globalTrack())];
          else if (muon->innerTrack().isNonnull())
            id = &leptons[reco::TrackBaseRef(muon->innerTrack())];
          else if (muon->outerTrack().isNonnull())
            // does this makes sense ?
            id = &leptons[reco::TrackBaseRef(muon->outerTrack())];
          else
            continue;
          if (haveLeptonCands)
            (*id)[SoftLeptonProperties::quality::btagMuonCands] = (*h_leptonCands)[h_muons->refAt(muon - h_muons->begin())];
        }
      }
      break;
    }
  } { // else
    // look for edm::View<Track> 
    Handle<edm::View<reco::Track> > h_tracks;
    event.getByLabel(m_leptons, h_tracks);
    if (h_tracks.isValid()) {
      for (unsigned int i = 0; i < h_tracks->size(); i++) {
        LeptonIds &id = leptons[h_tracks->refAt(i)];
        if (haveLeptonCands)
          id[SoftLeptonProperties::quality::btagLeptonCands] = (*h_leptonCands)[h_tracks->refAt(i)];
      }
      break;
    }
  } { // else
    throw edm::Exception(edm::errors::NotFound) << "Object " << m_leptons << " of type among (\"edm::View<reco::GsfElectron>\", \"edm::View<reco::Muon>\", \"edm::View<reco::Track>\") !found";
  } } while (false);

  if (!(m_leptonId == edm::InputTag())) {
    Handle< ValueMap<float> > h_leptonId;
    event.getByLabel(m_leptonId, h_leptonId);

    for (Leptons::iterator lepton = leptons.begin(); lepton != leptons.end(); ++lepton)
      lepton->second[leptonId] = (*h_leptonId)[lepton->first];
  }

  // output collections
  std::auto_ptr<reco::SoftLeptonTagInfoCollection> outputCollection(  new reco::SoftLeptonTagInfoCollection() );
  for (unsigned int i = 0; i < jets.size(); ++i) {
    reco::SoftLeptonTagInfo result = tag( jets[i], tracks[i], leptons, vertex );
    outputCollection->push_back( result );
  }
  event.put( outputCollection );
}
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 329 of file SoftLepton.cc.

References AXIS_CALORIMETRIC_NOLEPTON, 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_PI, m_refineJetAxis, reco::TrackBase::momentum(), perp(), reco::TrackBase::phi(), and reco::TrackBase::pt().

Referenced by tag().

        {
  math::XYZVector axis = jet->momentum();

  if (m_refineJetAxis == AXIS_CHARGED_AVERAGE or
      m_refineJetAxis == AXIS_CHARGED_AVERAGE_NOLEPTON) {

    double sum_pT        = 0.;
    double sum_eta_by_pT = 0.;
    double sum_phi_by_pT = 0.;

    double perp;
    double phi_rel;
    double eta_rel;

    // refine jet eta and phi with charged tracks measurements, if available
    for (reco::TrackRefVector::const_iterator track_it = tracks.begin(); track_it != tracks.end(); ++track_it ) {
      const reco::Track & track = **track_it;

      perp = track.pt();
      eta_rel = (double) track.eta() - axis.eta();
      phi_rel = (double) track.phi() - axis.phi();
      while (phi_rel < -M_PI) phi_rel += 2*M_PI;
      while (phi_rel >  M_PI) phi_rel -= 2*M_PI;

      sum_pT        += perp;
      sum_phi_by_pT += perp * phi_rel;
      sum_eta_by_pT += perp * eta_rel;
    }

    // "remove" excluded track
    if (m_refineJetAxis == AXIS_CHARGED_AVERAGE_NOLEPTON and exclude.isNonnull()) {
      const reco::Track & track = *exclude;

      perp = track.pt();
      eta_rel = (double) track.eta() - axis.eta();
      phi_rel = (double) track.phi() - axis.phi();
      while (phi_rel < -M_PI) phi_rel += 2*M_PI;
      while (phi_rel >  M_PI) phi_rel -= 2*M_PI;

      sum_pT        -= perp;
      sum_phi_by_pT -= perp * phi_rel;
      sum_eta_by_pT -= perp * eta_rel;
    }

    if (sum_pT > 1.)    // avoid the case of only the lepton-track with small rounding errors
      axis = math::RhoEtaPhiVector( axis.rho(), axis.eta() + sum_eta_by_pT / sum_pT, axis.phi() + sum_phi_by_pT / sum_pT);
    
  } else if (m_refineJetAxis == AXIS_CHARGED_SUM or
             m_refineJetAxis == AXIS_CHARGED_SUM_NOLEPTON) {
    math::XYZVector sum;

    // recalculate the jet direction as the sum of charget tracks momenta
    for (reco::TrackRefVector::const_iterator track_it = tracks.begin(); track_it != tracks.end(); ++track_it ) {
      const reco::Track & track = **track_it;
      sum += track.momentum();
    }

    // "remove" excluded track
    if (m_refineJetAxis == AXIS_CHARGED_SUM_NOLEPTON and exclude.isNonnull()) {
      const reco::Track & track = *exclude;
      sum -= track.momentum();
    }

    if (sum.R() > 1.) // avoid the case of only the lepton-track with small rounding errors
      axis = sum;
  } else if(m_refineJetAxis == AXIS_CALORIMETRIC_NOLEPTON) {
    axis -= exclude->momentum();
  }
  
  return GlobalVector(axis.x(), axis.y(), axis.z());
}
double SoftLepton::relativeEta ( const math::XYZVector vector,
const math::XYZVector axis 
) [static, protected]

Definition at line 405 of file SoftLepton.cc.

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

Referenced by tag().

                                                                                     {
  double mag = vector.r() * axis.r();
  double dot = vector.Dot(axis); 
  return -log((mag - dot)/(mag + dot)) / 2;
}
reco::SoftLeptonTagInfo SoftLepton::tag ( const edm::RefToBase< reco::Jet > &  jet,
const reco::TrackRefVector tracks,
const Leptons leptons,
const reco::Vertex primaryVertex 
) const

Definition at line 287 of file SoftLepton.cc.

References boostedPPar(), TransientTrackBuilder::build(), reco::SoftLeptonProperties::deltaR, reco::SoftLeptonProperties::etaRel, info, reco::SoftLeptonTagInfo::insert(), m_chi2Cut, m_deltaRCut, m_transientTrackBuilder, reco::SoftLeptonProperties::p0Par, reco::SoftLeptonProperties::ptRel, reco::SoftLeptonProperties::ratio, reco::SoftLeptonProperties::ratioRel, refineJetAxis(), relativeEta(), reco::JetTagInfo::setJetRef(), reco::SoftLeptonProperties::setQuality(), IPTools::signedImpactParameter3D(), IPTools::signedTransverseImpactParameter(), reco::SoftLeptonProperties::sip2d, reco::SoftLeptonProperties::sip3d, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by produce().

        {
  reco::SoftLeptonTagInfo info;
  info.setJetRef( jet );

  for(Leptons::const_iterator lepton = leptons.begin(); lepton != leptons.end(); ++lepton) {
    const math::XYZVector & lepton_momentum = lepton->first->momentum();
    if (m_chi2Cut > 0.0 && lepton->first->normalizedChi2() > m_chi2Cut)
      continue;

    const GlobalVector jetAxis = refineJetAxis( jet, tracks, lepton->first );
    const math::XYZVector axis( jetAxis.x(), jetAxis.y(), jetAxis.z());
    if (DeltaR(lepton_momentum, axis) > m_deltaRCut)
      continue;

    reco::SoftLeptonProperties properties;

    reco::TransientTrack transientTrack = m_transientTrackBuilder->build(*lepton->first);
    properties.sip2d    = IPTools::signedTransverseImpactParameter( transientTrack, jetAxis, primaryVertex ).second.significance();
    properties.sip3d    = IPTools::signedImpactParameter3D( transientTrack, jetAxis, primaryVertex ).second.significance();
    properties.deltaR   = DeltaR( lepton_momentum, axis );
    properties.ptRel    = Perp( lepton_momentum, axis );
    properties.p0Par    = boostedPPar( lepton_momentum, axis );
    properties.etaRel   = relativeEta( lepton_momentum, axis );
    properties.ratio    = lepton_momentum.R() / axis.R();
    properties.ratioRel = lepton_momentum.Dot(axis) / axis.Mag2();
 
    for(LeptonIds::const_iterator iter = lepton->second.begin(); iter != lepton->second.end(); ++iter)
      properties.setQuality(static_cast<SoftLeptonProperties::quality::Generic>(iter->first), iter->second);

    info.insert( lepton->first, properties );
  }

  return info;
}

Member Data Documentation

double SoftLepton::m_chi2Cut [private]

Definition at line 103 of file SoftLepton.h.

Referenced by tag().

double SoftLepton::m_deltaRCut [private]

Definition at line 102 of file SoftLepton.h.

Referenced by tag().

Definition at line 91 of file SoftLepton.h.

Referenced by produce().

Definition at line 94 of file SoftLepton.h.

Referenced by produce().

Definition at line 95 of file SoftLepton.h.

Referenced by produce().

Definition at line 93 of file SoftLepton.h.

Referenced by produce().

Definition at line 106 of file SoftLepton.h.

Referenced by produce().

Definition at line 92 of file SoftLepton.h.

Referenced by produce().

unsigned int SoftLepton::m_refineJetAxis [private]

Definition at line 101 of file SoftLepton.h.

Referenced by refineJetAxis().

Definition at line 98 of file SoftLepton.h.

Referenced by produce(), and tag().

Definition at line 109 of file SoftLepton.h.

Referenced by produce().