CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
EleTkIsolFromCands Class Reference

#include <EleTkIsolFromCands.h>

Classes

struct  TrkCuts
 

Public Member Functions

std::pair< int, double > calIsol (const reco::TrackBase &trk, const pat::PackedCandidateCollection &cands, const edm::View< reco::GsfElectron > &eles)
 
std::pair< int, double > calIsol (const double eleEta, const double elePhi, const double eleVZ, const pat::PackedCandidateCollection &cands, const edm::View< reco::GsfElectron > &eles)
 
double calIsolPt (const reco::TrackBase &trk, const pat::PackedCandidateCollection &cands, const edm::View< reco::GsfElectron > &eles)
 
double calIsolPt (const double eleEta, const double elePhi, const double eleVZ, const pat::PackedCandidateCollection &cands, const edm::View< reco::GsfElectron > &eles)
 
 EleTkIsolFromCands (const edm::ParameterSet &para)
 
 EleTkIsolFromCands (const EleTkIsolFromCands &)=default
 
EleTkIsolFromCandsoperator= (const EleTkIsolFromCands &)=default
 
 ~EleTkIsolFromCands ()=default
 

Static Public Member Functions

static bool passTrkSel (const reco::Track &trk, const double trkPt, const TrkCuts &cuts, const double eleEta, const double elePhi, const double eleVZ)
 
static edm::ParameterSetDescription pSetDescript ()
 

Private Member Functions

double getTrkPt (const reco::TrackBase &trk, const edm::View< reco::GsfElectron > &eles)
 

Static Private Member Functions

static bool passAlgo (const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackAlgorithm > &algosToRej)
 
static bool passQual (const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackQuality > &quals)
 

Private Attributes

TrkCuts barrelCuts_
 
TrkCuts endcapCuts_
 

Detailed Description

Definition at line 11 of file EleTkIsolFromCands.h.

Constructor & Destructor Documentation

EleTkIsolFromCands::EleTkIsolFromCands ( const edm::ParameterSet para)
explicit

Definition at line 47 of file EleTkIsolFromCands.cc.

47  :
48  barrelCuts_(para.getParameter<edm::ParameterSet>("barrelCuts")),
49  endcapCuts_(para.getParameter<edm::ParameterSet>("endcapCuts"))
50 {
51 
52 
53 }
T getParameter(std::string const &) const
EleTkIsolFromCands::EleTkIsolFromCands ( const EleTkIsolFromCands )
default
EleTkIsolFromCands::~EleTkIsolFromCands ( )
default

Member Function Documentation

std::pair< int, double > EleTkIsolFromCands::calIsol ( const reco::TrackBase trk,
const pat::PackedCandidateCollection cands,
const edm::View< reco::GsfElectron > &  eles 
)

Definition at line 64 of file EleTkIsolFromCands.cc.

References reco::TrackBase::eta(), reco::TrackBase::phi(), and reco::TrackBase::vz().

Referenced by calIsolPt().

67 {
68  return calIsol(eleTrk.eta(),eleTrk.phi(),eleTrk.vz(),cands,eles);
69 }
std::pair< int, double > calIsol(const reco::TrackBase &trk, const pat::PackedCandidateCollection &cands, const edm::View< reco::GsfElectron > &eles)
std::pair< int, double > EleTkIsolFromCands::calIsol ( const double  eleEta,
const double  elePhi,
const double  eleVZ,
const pat::PackedCandidateCollection cands,
const edm::View< reco::GsfElectron > &  eles 
)

Definition at line 72 of file EleTkIsolFromCands.cc.

References funct::abs(), barrelCuts_, particleFlowClusterECALTimeSelected_cfi::cuts, endcapCuts_, getTrkPt(), passTrkSel(), and reco::TrackBase::pt().

76 {
77 
78  double ptSum=0.;
79  int nrTrks=0;
80 
81  const TrkCuts& cuts = std::abs(eleEta)<1.5 ? barrelCuts_ : endcapCuts_;
82 
83  for(auto& cand : cands){
84  if(cand.charge()!=0){
85  const reco::Track& trk = cand.pseudoTrack();
86  double trkPt = std::abs(cand.pdgId())!=11 ? trk.pt() : getTrkPt(trk,eles);
87  if(passTrkSel(trk,trkPt,cuts,eleEta,elePhi,eleVZ)){
88  ptSum+=trkPt;
89  nrTrks++;
90  }
91  }
92  }
93  return {nrTrks,ptSum};
94 }
static bool passTrkSel(const reco::Track &trk, const double trkPt, const TrkCuts &cuts, const double eleEta, const double elePhi, const double eleVZ)
double pt() const
track transverse momentum
Definition: TrackBase.h:621
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double getTrkPt(const reco::TrackBase &trk, const edm::View< reco::GsfElectron > &eles)
double EleTkIsolFromCands::calIsolPt ( const reco::TrackBase trk,
const pat::PackedCandidateCollection cands,
const edm::View< reco::GsfElectron > &  eles 
)
inline

Definition at line 45 of file EleTkIsolFromCands.h.

References calIsol().

Referenced by ElectronHEEPIDValueMapProducer::calTrkIso().

46  {
47  return calIsol(trk,cands,eles).second;
48  }
std::pair< int, double > calIsol(const reco::TrackBase &trk, const pat::PackedCandidateCollection &cands, const edm::View< reco::GsfElectron > &eles)
double EleTkIsolFromCands::calIsolPt ( const double  eleEta,
const double  elePhi,
const double  eleVZ,
const pat::PackedCandidateCollection cands,
const edm::View< reco::GsfElectron > &  eles 
)
inline

Definition at line 50 of file EleTkIsolFromCands.h.

References calIsol(), particleFlowClusterECALTimeSelected_cfi::cuts, getTrkPt(), passAlgo(), passQual(), and passTrkSel().

52  {
53  return calIsol(eleEta,elePhi,eleVZ,cands,eles).second;
54  }
std::pair< int, double > calIsol(const reco::TrackBase &trk, const pat::PackedCandidateCollection &cands, const edm::View< reco::GsfElectron > &eles)
double EleTkIsolFromCands::getTrkPt ( const reco::TrackBase trk,
const edm::View< reco::GsfElectron > &  eles 
)
private

Definition at line 143 of file EleTkIsolFromCands.cc.

References funct::abs(), reco::deltaPhi(), reco::TrackBase::eta(), match(), reco::TrackBase::phi(), and reco::TrackBase::pt().

Referenced by calIsol(), calIsolPt(), and passAlgo().

145 {
146  //note, the trk.eta(),trk.phi() should be identical to the gsf track eta,phi
147  //although this may not be the case due to roundings after packing
148  auto match=[](const reco::TrackBase& trk,const reco::GsfElectron& ele){
149  return std::abs(trk.eta()-ele.gsfTrack()->eta())<0.001 &&
150  reco::deltaPhi(trk.phi(),ele.gsfTrack()->phi())<0.001;// &&
151  };
152  for(auto& ele : eles){
153  if(ele.gsfTrack().isNonnull()){
154  if(match(trk,ele)){
155  return ele.gsfTrack()->pt();
156  }
157  }
158  }
159  return trk.pt();
160 
161 }
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
double pt() const
track transverse momentum
Definition: TrackBase.h:621
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
EleTkIsolFromCands& EleTkIsolFromCands::operator= ( const EleTkIsolFromCands )
default
bool EleTkIsolFromCands::passAlgo ( const reco::TrackBase trk,
const std::vector< reco::TrackBase::TrackAlgorithm > &  algosToRej 
)
staticprivate

Definition at line 133 of file EleTkIsolFromCands.cc.

References reco::TrackBase::algo(), and getTrkPt().

Referenced by calIsolPt(), passQual(), and passTrkSel().

135 {
136  return algosToRej.empty() || !std::binary_search(algosToRej.begin(),algosToRej.end(),trk.algo());
137 }
TrackAlgorithm algo() const
Definition: TrackBase.h:497
bool EleTkIsolFromCands::passQual ( const reco::TrackBase trk,
const std::vector< reco::TrackBase::TrackQuality > &  quals 
)
staticprivate

Definition at line 120 of file EleTkIsolFromCands.cc.

References passAlgo(), and reco::TrackBase::quality().

Referenced by calIsolPt(), and passTrkSel().

122 {
123  if(quals.empty()) return true;
124 
125  for(auto qual : quals) {
126  if(trk.quality(qual)) return true;
127  }
128 
129  return false;
130 }
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:510
bool EleTkIsolFromCands::passTrkSel ( const reco::Track trk,
const double  trkPt,
const TrkCuts cuts,
const double  eleEta,
const double  elePhi,
const double  eleVZ 
)
static

Definition at line 97 of file EleTkIsolFromCands.cc.

References funct::abs(), EleTkIsolFromCands::TrkCuts::algosToReject, EleTkIsolFromCands::TrkCuts::allowedQualities, reco::deltaR2(), reco::TrackBase::eta(), reco::TrackBase::hitPattern(), EleTkIsolFromCands::TrkCuts::maxDPtPt, EleTkIsolFromCands::TrkCuts::maxDR2, EleTkIsolFromCands::TrkCuts::maxDZ, EleTkIsolFromCands::TrkCuts::minDEta, EleTkIsolFromCands::TrkCuts::minDR2, EleTkIsolFromCands::TrkCuts::minHits, EleTkIsolFromCands::TrkCuts::minPixelHits, EleTkIsolFromCands::TrkCuts::minPt, reco::HitPattern::numberOfValidHits(), reco::HitPattern::numberOfValidPixelHits(), passAlgo(), passQual(), reco::TrackBase::phi(), reco::TrackBase::ptError(), and reco::TrackBase::vz().

Referenced by calIsol(), and calIsolPt().

101 {
102  const float dR2 = reco::deltaR2(eleEta,elePhi,trk.eta(),trk.phi());
103  const float dEta = trk.eta()-eleEta;
104  const float dZ = eleVZ - trk.vz();
105 
106  return dR2>=cuts.minDR2 && dR2<=cuts.maxDR2 &&
107  std::abs(dEta)>=cuts.minDEta &&
108  std::abs(dZ)<cuts.maxDZ &&
109  trk.hitPattern().numberOfValidHits() >= cuts.minHits &&
110  trk.hitPattern().numberOfValidPixelHits() >=cuts.minPixelHits &&
111  (trk.ptError()/trkPt < cuts.maxDPtPt || cuts.maxDPtPt<0) &&
112  passQual(trk,cuts.allowedQualities) &&
113  passAlgo(trk,cuts.algosToReject) &&
114  trkPt > cuts.minPt;
115 }
int numberOfValidHits() const
Definition: HitPattern.h:823
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
static bool passAlgo(const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackAlgorithm > &algosToRej)
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:763
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static bool passQual(const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackQuality > &quals)
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:669
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:446
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
int numberOfValidPixelHits() const
Definition: HitPattern.h:838
edm::ParameterSetDescription EleTkIsolFromCands::pSetDescript ( )
static

Definition at line 55 of file EleTkIsolFromCands.cc.

References edm::ParameterSetDescription::add(), and EleTkIsolFromCands::TrkCuts::pSetDescript().

Referenced by ElectronHEEPIDValueMapProducer::fillDescriptions().

56 {
58  desc.add("barrelCuts",TrkCuts::pSetDescript());
59  desc.add("endcapCuts",TrkCuts::pSetDescript());
60  return desc;
61 }
static edm::ParameterSetDescription pSetDescript()
ParameterDescriptionBase * add(U const &iLabel, T const &value)

Member Data Documentation

TrkCuts EleTkIsolFromCands::barrelCuts_
private

Definition at line 29 of file EleTkIsolFromCands.h.

Referenced by calIsol().

TrkCuts EleTkIsolFromCands::endcapCuts_
private

Definition at line 29 of file EleTkIsolFromCands.h.

Referenced by calIsol().