CMS 3D CMS Logo

EleTkIsolFromCands.cc
Go to the documentation of this file.
4 
6  minPt = para.getParameter<double>("minPt");
7  auto sq = [](double val) { return val * val; };
8  minDR2 = sq(para.getParameter<double>("minDR"));
9  maxDR2 = sq(para.getParameter<double>("maxDR"));
10  minDEta = para.getParameter<double>("minDEta");
11  maxDZ = para.getParameter<double>("maxDZ");
12  minHits = para.getParameter<int>("minHits");
13  minPixelHits = para.getParameter<int>("minPixelHits");
14  maxDPtPt = para.getParameter<double>("maxDPtPt");
15 
16  auto qualNames = para.getParameter<std::vector<std::string> >("allowedQualities");
17  auto algoNames = para.getParameter<std::vector<std::string> >("algosToReject");
18 
19  for (auto& qualName : qualNames) {
21  }
22  for (auto& algoName : algoNames) {
24  }
25  std::sort(algosToReject.begin(), algosToReject.end());
26 }
27 
30  desc.add<double>("minPt", 1.0);
31  desc.add<double>("maxDR", 0.3);
32  desc.add<double>("minDR", 0.000);
33  desc.add<double>("minDEta", 0.005);
34  desc.add<double>("maxDZ", 0.1);
35  desc.add<double>("maxDPtPt", -1);
36  desc.add<int>("minHits", 8);
37  desc.add<int>("minPixelHits", 1);
38  desc.add<std::vector<std::string> >("allowedQualities");
39  desc.add<std::vector<std::string> >("algosToReject");
40  return desc;
41 }
42 
44  : barrelCuts_(para.getParameter<edm::ParameterSet>("barrelCuts")),
45  endcapCuts_(para.getParameter<edm::ParameterSet>("endcapCuts")) {}
46 
49  desc.add("barrelCuts", TrkCuts::pSetDescript());
50  desc.add("endcapCuts", TrkCuts::pSetDescript());
51  return desc;
52 }
53 
54 std::pair<int, double> EleTkIsolFromCands::calIsol(const reco::TrackBase& eleTrk,
56  const PIDVeto pidVeto) const {
57  return calIsol(eleTrk.eta(), eleTrk.phi(), eleTrk.vz(), cands, pidVeto);
58 }
59 
60 std::pair<int, double> EleTkIsolFromCands::calIsol(const double eleEta,
61  const double elePhi,
62  const double eleVZ,
64  const PIDVeto pidVeto) const {
65  double ptSum = 0.;
66  int nrTrks = 0;
67 
68  const TrkCuts& cuts = std::abs(eleEta) < 1.5 ? barrelCuts_ : endcapCuts_;
69 
70  for (auto& cand : cands) {
71  if (cand.hasTrackDetails() && cand.charge() != 0 && passPIDVeto(cand.pdgId(), pidVeto)) {
72  const reco::Track& trk = cand.pseudoTrack();
73  if (passTrkSel(trk, trk.pt(), cuts, eleEta, elePhi, eleVZ)) {
74  ptSum += trk.pt();
75  nrTrks++;
76  }
77  }
78  }
79  return {nrTrks, ptSum};
80 }
81 
82 std::pair<int, double> EleTkIsolFromCands::calIsol(const reco::TrackBase& eleTrk,
83  const reco::TrackCollection& tracks) const {
84  return calIsol(eleTrk.eta(), eleTrk.phi(), eleTrk.vz(), tracks);
85 }
86 
87 std::pair<int, double> EleTkIsolFromCands::calIsol(const double eleEta,
88  const double elePhi,
89  const double eleVZ,
90  const reco::TrackCollection& tracks) const {
91  double ptSum = 0.;
92  int nrTrks = 0;
93 
94  const TrkCuts& cuts = std::abs(eleEta) < 1.5 ? barrelCuts_ : endcapCuts_;
95 
96  for (auto& trk : tracks) {
97  if (passTrkSel(trk, trk.pt(), cuts, eleEta, elePhi, eleVZ)) {
98  ptSum += trk.pt();
99  nrTrks++;
100  }
101  }
102  return {nrTrks, ptSum};
103 }
104 
106  int pidAbs = std::abs(pdgId);
107  switch (veto) {
108  case PIDVeto::NONE:
109  return true;
110  case PIDVeto::ELES:
111  if (pidAbs == 11)
112  return false;
113  else
114  return true;
115  case PIDVeto::NONELES:
116  if (pidAbs == 11)
117  return true;
118  else
119  return false;
120  }
121  throw cms::Exception("CodeError") << "invalid PIDVeto " << static_cast<int>(veto) << ", "
122  << "this is likely due to some static casting of invalid ints somewhere";
123 }
124 
126  if (vetoStr == "NONE")
127  return PIDVeto::NONE;
128  else if (vetoStr == "ELES")
129  return PIDVeto::ELES;
130  else if (vetoStr == "NONELES")
131  return PIDVeto::NONELES;
132  else {
133  throw cms::Exception("CodeError") << "unrecognised string " << vetoStr
134  << ", either a typo or this function needs to be updated";
135  }
136 }
137 
139  const double trkPt,
140  const TrkCuts& cuts,
141  const double eleEta,
142  const double elePhi,
143  const double eleVZ) {
144  const float dR2 = reco::deltaR2(eleEta, elePhi, trk.eta(), trk.phi());
145  const float dEta = trk.eta() - eleEta;
146  const float dZ = eleVZ - trk.vz();
147 
148  return dR2 >= cuts.minDR2 && dR2 <= cuts.maxDR2 && std::abs(dEta) >= cuts.minDEta && std::abs(dZ) < cuts.maxDZ &&
149  trk.hitPattern().numberOfValidHits() >= cuts.minHits &&
150  trk.hitPattern().numberOfValidPixelHits() >= cuts.minPixelHits &&
151  (trk.ptError() / trkPt < cuts.maxDPtPt || cuts.maxDPtPt < 0) && passQual(trk, cuts.allowedQualities) &&
152  passAlgo(trk, cuts.algosToReject) && trkPt > cuts.minPt;
153 }
154 
155 bool EleTkIsolFromCands::passQual(const reco::TrackBase& trk, const std::vector<reco::TrackBase::TrackQuality>& quals) {
156  if (quals.empty())
157  return true;
158 
159  for (auto qual : quals) {
160  if (trk.quality(qual))
161  return true;
162  }
163 
164  return false;
165 }
166 
168  const std::vector<reco::TrackBase::TrackAlgorithm>& algosToRej) {
169  return algosToRej.empty() || !std::binary_search(algosToRej.begin(), algosToRej.end(), trk.algo());
170 }
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
EleTkIsolFromCands::TrkCuts::maxDPtPt
float maxDPtPt
Definition: EleTkIsolFromCands.h:61
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
reco::TrackBase::ptError
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:702
EleTkIsolFromCands::passPIDVeto
static bool passPIDVeto(const int pdgId, const EleTkIsolFromCands::PIDVeto pidVeto)
Definition: EleTkIsolFromCands.cc:105
EleTkIsolFromCands::TrkCuts::minHits
float minHits
Definition: EleTkIsolFromCands.h:59
EleTkIsolFromCands::barrelCuts_
TrkCuts barrelCuts_
Definition: EleTkIsolFromCands.h:68
EleTkIsolFromCands::TrkCuts::allowedQualities
std::vector< reco::TrackBase::TrackQuality > allowedQualities
Definition: EleTkIsolFromCands.h:62
edm
HLT enums.
Definition: AlignableModifier.h:19
EleTkIsolFromCands::PIDVeto::NONE
EleTkIsolFromCands::pidVetoFromStr
static PIDVeto pidVetoFromStr(const std::string &vetoStr)
Definition: EleTkIsolFromCands.cc:125
EleTkIsolFromCands::TrkCuts::pSetDescript
static edm::ParameterSetDescription pSetDescript()
Definition: EleTkIsolFromCands.cc:28
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
EleTkIsolFromCands::PIDVeto::NONELES
HIPAlignmentAlgorithm_cfi.algoName
algoName
Definition: HIPAlignmentAlgorithm_cfi.py:5
EleTkIsolFromCands::calIsol
std::pair< int, double > calIsol(const reco::TrackBase &trk, const pat::PackedCandidateCollection &cands, const PIDVeto=PIDVeto::NONE) const
Definition: EleTkIsolFromCands.cc:54
EleTkIsolFromCands::PIDVeto::ELES
HLT_2018_cff.dEta
dEta
Definition: HLT_2018_cff.py:12289
deltaR.h
reco::TrackBase::pt
double pt() const
track transverse momentum
Definition: TrackBase.h:608
EleTkIsolFromCands::TrkCuts::TrkCuts
TrkCuts(const edm::ParameterSet &para)
Definition: EleTkIsolFromCands.cc:5
EleTkIsolFromCands::PIDVeto
PIDVeto
Definition: EleTkIsolFromCands.h:46
Track.h
EleTkIsolFromCands::TrkCuts::maxDR2
float maxDR2
Definition: EleTkIsolFromCands.h:56
reco::TrackBase::vz
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:632
reco::Track
Definition: Track.h:27
reco::TrackBase::phi
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:620
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
EleTkIsolFromCands::passQual
static bool passQual(const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackQuality > &quals)
Definition: EleTkIsolFromCands.cc:155
edm::ParameterSet
Definition: ParameterSet.h:36
EleTkIsolFromCands::TrkCuts::minDR2
float minDR2
Definition: EleTkIsolFromCands.h:55
EleTkIsolFromCands.h
ParameterSet
Definition: Functions.h:16
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
reco::TrackBase::eta
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:623
cand
Definition: decayParser.h:34
reco::TrackBase
Definition: TrackBase.h:62
EleTkIsolFromCands::passAlgo
static bool passAlgo(const reco::TrackBase &trk, const std::vector< reco::TrackBase::TrackAlgorithm > &algosToRej)
Definition: EleTkIsolFromCands.cc:167
EleTkIsolFromCands::pSetDescript
static edm::ParameterSetDescription pSetDescript()
Definition: EleTkIsolFromCands.cc:47
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:118
pat::PackedCandidateCollection
std::vector< pat::PackedCandidate > PackedCandidateCollection
Definition: PackedCandidate.h:1130
reco::TrackBase::algo
TrackAlgorithm algo() const
Definition: TrackBase.h:532
EleTkIsolFromCands::TrkCuts::minPt
float minPt
Definition: EleTkIsolFromCands.h:54
reco::TrackBase::qualityByName
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
reco::TrackBase::hitPattern
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:489
heppy_batch.val
val
Definition: heppy_batch.py:351
EleTkIsolFromCands::passTrkSel
static bool passTrkSel(const reco::TrackBase &trk, const double trkPt, const TrkCuts &cuts, const double eleEta, const double elePhi, const double eleVZ)
Definition: EleTkIsolFromCands.cc:138
Exception
Definition: hltDiff.cc:246
EleTkIsolFromCands::TrkCuts
Definition: EleTkIsolFromCands.h:53
EleTkIsolFromCands::TrkCuts::algosToReject
std::vector< reco::TrackBase::TrackAlgorithm > algosToReject
Definition: EleTkIsolFromCands.h:63
HLT_2018_cff.cands
cands
Definition: HLT_2018_cff.py:13762
EleTkIsolFromCands::EleTkIsolFromCands
EleTkIsolFromCands(const edm::ParameterSet &para)
Definition: EleTkIsolFromCands.cc:43
EleTkIsolFromCands::TrkCuts::maxDZ
float maxDZ
Definition: EleTkIsolFromCands.h:58
reco::HitPattern::numberOfValidPixelHits
int numberOfValidPixelHits() const
Definition: HitPattern.h:800
reco::TrackBase::algoByName
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:137
L1TMuonDQMOffline_cfi.cuts
cuts
Definition: L1TMuonDQMOffline_cfi.py:41
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::HitPattern::numberOfValidHits
int numberOfValidHits() const
Definition: HitPattern.h:786
PbPb_ZMuSkimMuonDPG_cff.veto
veto
Definition: PbPb_ZMuSkimMuonDPG_cff.py:61
EleTkIsolFromCands::endcapCuts_
TrkCuts endcapCuts_
Definition: EleTkIsolFromCands.h:68
reco::TrackBase::quality
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:537
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
EleTkIsolFromCands::TrkCuts::minDEta
float minDEta
Definition: EleTkIsolFromCands.h:57
EleTkIsolFromCands::TrkCuts::minPixelHits
float minPixelHits
Definition: EleTkIsolFromCands.h:60