CMS 3D CMS Logo

PATTauHybridProducer.cc
Go to the documentation of this file.
13 
15 public:
16  explicit PATTauHybridProducer(const edm::ParameterSet&);
17  ~PATTauHybridProducer() override{};
18 
19  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
20  void produce(edm::Event&, const edm::EventSetup&) override;
21 
22 private:
23  void fillTauFromJet(reco::PFTau& pfTau, const reco::JetBaseRef& jet);
24  //--- configuration parameters
29  const float dR2Max_, jetPtMin_, jetEtaMax_;
31  std::vector<std::string> pnetTauScoreNames_;
32  std::vector<std::string> pnetJetScoreNames_;
33  std::vector<std::string> pnetLepScoreNames_;
38 
39  const std::map<std::string, int> tagToDM_;
40  enum class tauId_pnet_idx : size_t { dm = 0, vsjet, vse, vsmu, ptcorr, qconf, pdm0, pdm1, pdm2, pdm10, pdm11, last };
41  enum class tauId_min_idx : size_t { hpsnew = 0, last };
42 };
43 
45  : tausToken_(consumes<pat::TauCollection>(cfg.getParameter<edm::InputTag>("src"))),
46  jetsToken_(consumes<pat::JetCollection>(cfg.getParameter<edm::InputTag>("jetSource"))),
47  addGenJetMatch_(cfg.getParameter<bool>("addGenJetMatch")),
48  dR2Max_(std::pow(cfg.getParameter<double>("dRMax"), 2)),
49  jetPtMin_(cfg.getParameter<double>("jetPtMin")),
50  jetEtaMax_(cfg.getParameter<double>("jetEtaMax")),
51  pnetLabel_(cfg.getParameter<std::string>("pnetLabel")),
52  pnetPtCorrName_(cfg.getParameter<std::string>("pnetPtCorrName")),
53  tauScoreMin_(cfg.getParameter<double>("tauScoreMin")),
54  vsJetMin_(cfg.getParameter<double>("vsJetMin")),
55  chargeAssignmentProbMin_(cfg.getParameter<double>("chargeAssignmentProbMin")),
56  checkTauScoreIsBest_(cfg.getParameter<bool>("checkTauScoreIsBest")),
57  usePFLeptonsAsChargedHadrons_(cfg.getParameter<bool>("usePFLeptonsAsChargedHadrons")),
58  tagToDM_({{"1h0p", 0}, {"1h1or2p", 1}, {"1h1p", 1}, {"1h2p", 2}, {"3h0p", 10}, {"3h1p", 11}}) {
59  // Read the different PNet score names
60  std::vector<std::string> pnetScoreNames = cfg.getParameter<std::vector<std::string>>("pnetScoreNames");
61  for (const auto& scoreName : pnetScoreNames) {
62  size_t labelLenght = scoreName.find(':') == std::string::npos ? 0 : scoreName.find(':') + 1;
63  std::string name = scoreName.substr(labelLenght);
64  if (name.find("prob") == std::string::npos)
65  continue;
66  if (name.find("tau") != std::string::npos)
67  pnetTauScoreNames_.push_back(name);
68  else if (name.find("ele") != std::string::npos || name.find("mu") != std::string::npos)
69  pnetLepScoreNames_.push_back(name);
70  else
71  pnetJetScoreNames_.push_back(name);
72  if (pnetPtCorrName_.find(':') != std::string::npos)
73  pnetPtCorrName_ = pnetPtCorrName_.substr(pnetPtCorrName_.find(':') + 1);
74  // GenJet matching
75  if (addGenJetMatch_) {
76  genJetMatchToken_ =
77  consumes<edm::Association<reco::GenJetCollection>>(cfg.getParameter<edm::InputTag>("genJetMatch"));
78  }
79  }
80 
81  produces<std::vector<pat::Tau>>();
82  //FIXME: produce a separate collection for PNet-recovered taus?
83 }
84 
86  // Get the vector of taus
88  evt.getByToken(tausToken_, inputTaus);
89 
90  auto outputTaus = std::make_unique<std::vector<pat::Tau>>();
91  outputTaus->reserve(inputTaus->size());
92 
93  // Get the vector of jets
96 
97  // Switch off gen-matching for real data
98  if (evt.isRealData()) {
99  addGenJetMatch_ = false;
100  }
102  if (addGenJetMatch_)
104 
105  // Minimal HPS-like tauID list
106  std::vector<pat::Tau::IdPair> tauIds_minimal((size_t)tauId_min_idx::last);
107  tauIds_minimal[(size_t)tauId_min_idx::hpsnew] = std::make_pair("decayModeFindingNewDMs", -1);
108 
109  // PNet tauID list
110  std::vector<pat::Tau::IdPair> tauIds_pnet((size_t)tauId_pnet_idx::last);
111  tauIds_pnet[(size_t)tauId_pnet_idx::dm] = std::make_pair("byPNetDecayMode", reco::PFTau::kNull);
112  tauIds_pnet[(size_t)tauId_pnet_idx::vsjet] = std::make_pair("byPNetVSjetraw", -1);
113  tauIds_pnet[(size_t)tauId_pnet_idx::vse] = std::make_pair("byPNetVSeraw", -1);
114  tauIds_pnet[(size_t)tauId_pnet_idx::vsmu] = std::make_pair("byPNetVSmuraw", -1);
115  tauIds_pnet[(size_t)tauId_pnet_idx::ptcorr] = std::make_pair("byPNetPtCorr", 1);
116  tauIds_pnet[(size_t)tauId_pnet_idx::qconf] = std::make_pair("byPNetQConf", 0);
117  tauIds_pnet[(size_t)tauId_pnet_idx::pdm0] = std::make_pair("byPNetProb1h0pi0", -1);
118  tauIds_pnet[(size_t)tauId_pnet_idx::pdm1] = std::make_pair("byPNetProb1h1pi0", -1);
119  tauIds_pnet[(size_t)tauId_pnet_idx::pdm2] = std::make_pair("byPNetProb1h2pi0", -1);
120  tauIds_pnet[(size_t)tauId_pnet_idx::pdm10] = std::make_pair("byPNetProb3h0pi0", -1);
121  tauIds_pnet[(size_t)tauId_pnet_idx::pdm11] = std::make_pair("byPNetProb3h1pi0", -1);
122 
123  std::set<unsigned int> matched_taus;
124  size_t jet_idx = 0;
125  for (const auto& jet : *jets) {
126  jet_idx++;
127  if (jet.pt() < jetPtMin_)
128  continue;
129  if (std::abs(jet.eta()) > jetEtaMax_)
130  continue;
131  size_t tau_idx = 0;
132  bool matched = false;
133 
134  // Analyse PNet scores
135  std::pair<std::string, float> bestPnetTauScore("probtauundef", -1);
136  float sumOfPnetTauScores = 0;
137  std::vector<float> tauPerDMScores(5);
138  float plusChargeProb = 0;
139  for (const auto& scoreName : pnetTauScoreNames_) {
140  float score = jet.bDiscriminator(pnetLabel_ + ":" + scoreName);
141  sumOfPnetTauScores += score;
142  if (scoreName.find("taup") != std::string::npos)
143  plusChargeProb += score;
144  if (score > bestPnetTauScore.second) {
145  bestPnetTauScore.second = score;
146  bestPnetTauScore.first = scoreName;
147  }
148  if (scoreName.find("1h0p") != std::string::npos)
149  tauPerDMScores[0] += score;
150  else if (scoreName.find("1h1") !=
151  std::string::
152  npos) //Note: final "p" in "1p" ommited to enble matching also with "1h1or2p" from early trainings
153  tauPerDMScores[1] += score;
154  else if (scoreName.find("1h2p") != std::string::npos)
155  tauPerDMScores[2] += score;
156  else if (scoreName.find("3h0p") != std::string::npos)
157  tauPerDMScores[3] += score;
158  else if (scoreName.find("3h1p") != std::string::npos)
159  tauPerDMScores[4] += score;
160  }
161  if (sumOfPnetTauScores > 0)
162  plusChargeProb /= sumOfPnetTauScores;
163 
164  float sumOfPnetEleScores = 0, sumOfPnetMuScores = 0;
165  bool isTauScoreBest = (sumOfPnetTauScores > 0);
166  for (const auto& scoreName : pnetLepScoreNames_) {
167  float score = jet.bDiscriminator(pnetLabel_ + ":" + scoreName);
168  if (scoreName.find("ele") != std::string::npos)
169  sumOfPnetEleScores += score;
170  else if (scoreName.find("mu") != std::string::npos)
171  sumOfPnetMuScores += score;
172  if (score > bestPnetTauScore.second)
173  isTauScoreBest = false;
174  }
175  if (checkTauScoreIsBest_ && isTauScoreBest) { //if needed iterate over jet scores
176  for (const auto& scoreName : pnetJetScoreNames_)
177  if (jet.bDiscriminator(pnetLabel_ + ":" + scoreName) > bestPnetTauScore.second)
178  isTauScoreBest = false;
179  }
180 
181  // PNet discriminants vs jets, electrons and muons
182  tauIds_pnet[(size_t)tauId_pnet_idx::vsjet].second =
183  sumOfPnetTauScores /
184  (1. - sumOfPnetEleScores -
185  sumOfPnetMuScores); //vsJet: tau scores by sum of tau and jet scores or equally by 1 - sum of lepton scores
186  tauIds_pnet[(size_t)tauId_pnet_idx::vse].second =
187  sumOfPnetTauScores /
188  (sumOfPnetTauScores + sumOfPnetEleScores); //vsEle: tau scores by sum of tau and ele scores
189  tauIds_pnet[(size_t)tauId_pnet_idx::vsmu].second =
190  sumOfPnetTauScores / (sumOfPnetTauScores + sumOfPnetMuScores); //vsMu: tau scores by sum of tau and mu scores
191 
192  // Decay mode and charge of the highest tau score
193  int bestCharge = 0;
194  size_t pos =
195  bestPnetTauScore.first.find("tau") + 3; //this is well defined by constraction as name is "probtauXXXX"
196  const char q = (pos < bestPnetTauScore.first.size()) ? bestPnetTauScore.first[pos] : 'u';
197  if (q == 'm') { //minus
198  pos++;
199  bestCharge = -1;
200  } else if (q == 'p') { //plus
201  pos++;
202  bestCharge = 1;
203  }
204  auto pNetDM = tagToDM_.find(bestPnetTauScore.first.substr(pos));
205  if (pNetDM != tagToDM_.end())
206  tauIds_pnet[(size_t)tauId_pnet_idx::dm].second = pNetDM->second;
207 
208  // PNet Pt correction
209  float ptcorr = jet.bDiscriminator(pnetLabel_ + ":" + pnetPtCorrName_);
210  if (ptcorr > -1000.) // -1000. is default for not found discriminantor
211  tauIds_pnet[(size_t)tauId_pnet_idx::ptcorr].second = ptcorr;
212 
213  // PNet charge confidence
214  tauIds_pnet[(size_t)tauId_pnet_idx::qconf].second = (plusChargeProb - 0.5);
215 
216  // PNet per decay mode normalised score
217  tauIds_pnet[(size_t)tauId_pnet_idx::pdm0].second = tauPerDMScores[0] / sumOfPnetTauScores;
218  tauIds_pnet[(size_t)tauId_pnet_idx::pdm1].second = tauPerDMScores[1] / sumOfPnetTauScores;
219  tauIds_pnet[(size_t)tauId_pnet_idx::pdm2].second = tauPerDMScores[2] / sumOfPnetTauScores;
220  tauIds_pnet[(size_t)tauId_pnet_idx::pdm10].second = tauPerDMScores[3] / sumOfPnetTauScores;
221  tauIds_pnet[(size_t)tauId_pnet_idx::pdm11].second = tauPerDMScores[4] / sumOfPnetTauScores;
222 
223  // Search for matching tau
224  for (const auto& inputTau : *inputTaus) {
225  tau_idx++;
226  if (matched_taus.count(tau_idx - 1) > 0)
227  continue;
228  float dR2 = deltaR2(jet, inputTau);
229  // select 1st found match rather than best match (both should be equivalent for reasonable dRMax)
230  if (dR2 < dR2Max_) {
231  matched_taus.insert(tau_idx - 1);
232  pat::Tau outputTau(inputTau);
233  const size_t nTauIds = inputTau.tauIDs().size();
234  std::vector<pat::Tau::IdPair> tauIds(nTauIds + tauIds_pnet.size());
235  for (size_t i = 0; i < nTauIds; ++i)
236  tauIds[i] = inputTau.tauIDs()[i];
237  for (size_t i = 0; i < tauIds_pnet.size(); ++i)
238  tauIds[nTauIds + i] = tauIds_pnet[i];
239  outputTau.setTauIDs(tauIds);
240  matched = true;
241  outputTaus->push_back(outputTau);
242 
243  break;
244  }
245  } // end of tau loop
246  if (matched)
247  continue;
248 
249  // Accept only jets passing minimal tau-like selection, i.e. with one of the tau score being globally the best and above some threshold, and with good quality of charge assignment
250  if ((checkTauScoreIsBest_ && !isTauScoreBest) || bestPnetTauScore.second < tauScoreMin_ ||
251  tauIds_pnet[(size_t)tauId_pnet_idx::vsjet].second < vsJetMin_ ||
252  std::abs(0.5 - plusChargeProb) < chargeAssignmentProbMin_)
253  continue;
254 
255  // Build taus from non-matched jets
256  // "Null" pftau with raw (uncorrected) jet kinematics
257  reco::PFTau pfTauFromJet(bestCharge, jet.correctedP4("Uncorrected"));
258  // Set PDGid
259  pfTauFromJet.setPdgId(bestCharge < 0 ? 15 : -15);
260  // and decay mode predicted by PNet
261  pfTauFromJet.setDecayMode(
262  static_cast<const reco::PFTau::hadronicDecayMode>(int(tauIds_pnet[(size_t)tauId_pnet_idx::dm].second)));
263  // Fill tau content using only jet consistunets within cone around leading
264  // charged particle
265  // FIXME: more sophisticated finding of tau constituents will be considered later
266  pfTauFromJet.setSignalConeSize(
267  std::clamp(3.6 / jet.correctedP4("Uncorrected").pt(), 0.08, 0.12)); // shrinking cone in function of jet-Pt
268  const edm::Ref<pat::JetCollection> jetRef(jets, jet_idx - 1);
269  fillTauFromJet(pfTauFromJet, reco::JetBaseRef(jetRef));
270 
271  // PATTau
272  pat::Tau outputTauFromJet(pfTauFromJet);
273  // Add tauIDs
274  std::vector<pat::Tau::IdPair> newtauIds(tauIds_minimal.size() + tauIds_pnet.size());
275  for (size_t i = 0; i < tauIds_minimal.size(); ++i)
276  newtauIds[i] = tauIds_minimal[i];
277  for (size_t i = 0; i < tauIds_pnet.size(); ++i)
278  newtauIds[tauIds_minimal.size() + i] = tauIds_pnet[i];
279  outputTauFromJet.setTauIDs(newtauIds);
280  // Add genTauJet match
281  if (addGenJetMatch_) {
282  reco::GenJetRef genJetTau = (*genJetMatch)[jetRef];
283  if (genJetTau.isNonnull() && genJetTau.isAvailable()) {
284  outputTauFromJet.setGenJet(genJetTau);
285  }
286  }
287  outputTaus->push_back(outputTauFromJet);
288 
289  } // end of jet loop
290 
291  // Taus non-matched to jets (usually at pt-threshold or/and eta boundaries)
292  if (matched_taus.size() < inputTaus->size()) {
293  for (size_t iTau = 0; iTau < inputTaus->size(); ++iTau) {
294  if (matched_taus.count(iTau) > 0)
295  continue;
296  const pat::Tau& inputTau = inputTaus->at(iTau);
297  pat::Tau outputTau(inputTau);
298  const size_t nTauIds = inputTau.tauIDs().size();
299  std::vector<pat::Tau::IdPair> tauIds(nTauIds + tauIds_pnet.size());
300  for (size_t i = 0; i < nTauIds; ++i)
301  tauIds[i] = inputTau.tauIDs()[i];
302  for (size_t i = 0; i < tauIds_pnet.size(); ++i) {
303  tauIds[nTauIds + i] = tauIds_pnet[i];
304  tauIds[nTauIds + i].second =
305  (i != (size_t)tauId_pnet_idx::ptcorr ? (i != (size_t)tauId_pnet_idx::qconf ? -1 : 0) : 1);
306  }
307  outputTau.setTauIDs(tauIds);
308  outputTaus->push_back(outputTau);
309  }
310  } //non-matched taus
311 
312  evt.put(std::move(outputTaus));
313 }
314 
316  // Use tools as in PFTau builders to select tau decay products and isolation candidates
317  typedef std::vector<reco::CandidatePtr> CandPtrs;
318 
319  // Get the charged hadron candidates
320  CandPtrs pfChs, pfChsSig;
321  // Check if we want to include electrons and muons in "charged hadron"
322  // collection. This is the preferred behavior, as the PF lepton selections
323  // are very loose.
325  pfChs = reco::tau::pfCandidatesByPdgId(*jet, 211);
326  } else {
327  pfChs = reco::tau::pfChargedCands(*jet);
328  }
329  // take 1st charged candidate with charge as of tau (collection is pt-sorted)
330  if (pfTau.charge() == 0 || pfChs.size() == 1) {
331  pfTau.setleadChargedHadrCand(pfChs[0]);
332  pfTau.setleadCand(pfChs[0]);
333  pfChsSig.push_back(pfChs[0]);
334  pfChs.erase(pfChs.begin());
335  } else {
336  for (CandPtrs::iterator it = pfChs.begin(); it != pfChs.end();) {
337  if ((*it)->charge() == pfTau.charge()) {
338  pfTau.setleadChargedHadrCand(*it);
339  pfTau.setleadCand(*it);
340  pfChsSig.push_back(*it);
341  it = pfChs.erase(it);
342  break;
343  } else {
344  ++it;
345  }
346  }
347  // In case of lack of candidate with charge same as of tau use leading charged candidate
348  if (pfTau.leadChargedHadrCand().isNull() && !pfChs.empty()) {
349  pfTau.setleadChargedHadrCand(pfChs[0]);
350  pfTau.setleadCand(pfChs[0]);
351  pfChsSig.push_back(pfChs[0]);
352  pfChs.erase(pfChs.begin());
353  }
354  }
355  // if more than one charged decay product is expected add all inside signal
356  // cone around the leading track
357  const double dR2Max = std::pow(pfTau.signalConeSize(), 2);
359  for (CandPtrs::iterator it = pfChs.begin(); it != pfChs.end();) {
360  if (deltaR2((*it)->p4(), pfTau.leadChargedHadrCand()->p4()) < dR2Max) {
361  pfChsSig.push_back(*it);
362  it = pfChs.erase(it);
363  } else {
364  ++it;
365  }
366  }
367  }
368  // Clean isolation candidates from low-pt and leptonic ones
369  for (CandPtrs::iterator it = pfChs.begin(); it != pfChs.end();) {
370  if ((*it)->pt() < 0.5 || std::abs((*it)->pdgId()) != 211) {
371  it = pfChs.erase(it);
372  } else {
373  ++it;
374  }
375  }
376  // Set charged candidates
377  pfTau.setsignalChargedHadrCands(pfChsSig);
378  pfTau.setisolationChargedHadrCands(pfChs);
379 
380  // Get the gamma candidates (pi0 decay products)
381  CandPtrs pfGammas, pfGammasSig;
383  // In case of lack of leading charged candidate substiute it with leading gamma candidate
384  if (pfTau.leadChargedHadrCand().isNull() && !pfGammas.empty()) {
386  pfTau.setleadCand(pfGammas[0]);
387  pfGammasSig.push_back(pfGammas[0]);
388  pfChs.erase(pfGammas.begin());
389  }
390  // Clean gamma candidates from low-pt ones
391  for (CandPtrs::iterator it = pfGammas.begin(); it != pfGammas.end();) {
392  if ((*it)->pt() < 0.5) {
393  it = pfGammas.erase(it);
394  } else {
395  ++it;
396  }
397  }
398  // if decay mode with pi0s is expected look for signal gamma candidates
399  // within eta-phi strips around leading track
400  if (pfTau.decayMode() % 5 != 0 && pfTau.leadChargedHadrCand().isNonnull()) {
401  for (CandPtrs::iterator it = pfGammas.begin(); it != pfGammas.end();) {
402  if (std::abs((*it)->eta() - pfTau.leadChargedHadrCand()->eta()) <
403  std::clamp(0.2 * std::pow((*it)->pt(), -0.66), 0.05, 0.15) &&
404  deltaPhi((*it)->phi(), pfTau.leadChargedHadrCand()->phi()) <
405  std::clamp(0.35 * std::pow((*it)->pt(), -0.71), 0.05, 0.3)) {
406  pfGammasSig.push_back(*it);
407  it = pfGammas.erase(it);
408  } else {
409  ++it;
410  }
411  }
412  }
413  // Set gamma candidates
414  pfTau.setsignalGammaCands(pfGammasSig);
416 }
417 
419  // patTauHybridProducer
421 
422  desc.add<edm::InputTag>("src", edm::InputTag("slimmedTaus"));
423  desc.add<edm::InputTag>("jetSource", edm::InputTag("slimmedJetsUpadted"));
424  desc.add<double>("dRMax", 0.4);
425  desc.add<double>("jetPtMin", 20.0);
426  desc.add<double>("jetEtaMax", 2.5);
427  desc.add<std::string>("pnetLabel", "pfParticleNetAK4JetTags");
428  desc.add<std::vector<std::string>>(
429  "pnetScoreNames",
430  {"probmu", "probele", "probtaup1h0p", "probtaup1h1p", "probtaup1h2p", "probtaup3h0p", "probtaup3h1p",
431  "probtaum1h0p", "probtaum1h1p", "probtaum1h2p", "probtaum3h0p", "probtaum3h1p", "probb", "probc",
432  "probuds", "probg", "ptcorr", "ptreshigh", "ptreslow", "ptnu"});
433  desc.add<std::string>("pnetPtCorrName", "ptcorr");
434  desc.add<double>("tauScoreMin", -1)->setComment("Minimal value of the best tau score to built recovery tau");
435  desc.add<double>("vsJetMin", -1)->setComment("Minimal value of PNet tau-vs-jet discriminant to built recovery tau");
436  desc.add<bool>("checkTauScoreIsBest", false)
437  ->setComment("If true, recovery tau is built only if one of tau scores is the highest");
438  desc.add<double>("chargeAssignmentProbMin", 0.2)
439  ->setComment("Minimal value of charge assignment probability to built recovery tau (0,0.5)");
440  desc.add<bool>("addGenJetMatch", true)->setComment("add MC genTauJet matching");
441  desc.add<edm::InputTag>("genJetMatch", edm::InputTag("tauGenJetMatch"));
442  desc.add<bool>("usePFLeptonsAsChargedHadrons", true)
443  ->setComment("If true, all charged particles are used as charged hadron candidates");
444 
445  descriptions.addWithDefaultLabel(desc);
446 }
447 
449 
std::vector< CandidatePtr > pfCandidatesByPdgId(const Jet &jet, int pdgId, bool sort=true)
const double dR2Max
void setsignalChargedHadrCands(const std::vector< reco::CandidatePtr > &)
Definition: PFTau.cc:77
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
double signalConeSize() const
Size of signal cone.
Definition: PFTau.h:174
std::vector< Jet > JetCollection
Definition: Jet.h:53
std::vector< std::string > pnetTauScoreNames_
const CandidatePtr & leadChargedHadrCand() const
Definition: PFTau.cc:63
edm::EDGetTokenT< pat::JetCollection > jetsToken_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
std::vector< Tau > TauCollection
Definition: Tau.h:35
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:540
void setTauIDs(const std::vector< IdPair > &ids)
Definition: Tau.h:352
void setleadCand(const CandidatePtr &)
Definition: PFTau.cc:69
std::vector< CandidatePtr > pfChargedCands(const Jet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
void setisolationChargedHadrCands(const std::vector< reco::CandidatePtr > &)
Definition: PFTau.cc:92
PATTauHybridProducer(const edm::ParameterSet &)
Definition: HeavyIon.h:7
const std::vector< IdPair > & tauIDs() const
Definition: Tau.h:349
void setleadChargedHadrCand(const CandidatePtr &)
Definition: PFTau.cc:67
U second(std::pair< T, U > const &p)
bool isNull() const
Checks for null.
Definition: Ptr.h:142
std::vector< reco::CandidatePtr > CandPtrs
std::vector< std::string > pnetJetScoreNames_
const float chargeAssignmentProbMin_
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
bool isAvailable() const
Definition: Ref.h:541
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const std::map< std::string, int > tagToDM_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
genJetMatch
switch on/off embedding of matched genJet&#39;s
Analysis-level tau class.
Definition: Tau.h:53
edm::EDGetTokenT< edm::Association< reco::GenJetCollection > > genJetMatchToken_
void produce(edm::Event &, const edm::EventSetup &) override
void fillTauFromJet(reco::PFTau &pfTau, const reco::JetBaseRef &jet)
void setGenJet(const reco::GenJetRef &ref)
set the matched GenJet
edm::EDGetTokenT< pat::TauCollection > tausToken_
const bool usePFLeptonsAsChargedHadrons_
hadronicDecayMode decayMode() const
Definition: PFTau.cc:325
void setisolationGammaCands(const std::vector< reco::CandidatePtr > &)
Definition: PFTau.cc:100
void setsignalGammaCands(const std::vector< reco::CandidatePtr > &)
Definition: PFTau.cc:85
HLT enums.
bool isRealData() const
Definition: EventBase.h:66
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< std::string > pnetLepScoreNames_
const std::string pnetLabel_
void setPdgId(int pdgId) final
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511
int charge() const final
electric charge
std::vector< CandidatePtr > pfGammas(const Jet &jet, bool sort=true)
Extract all pfGammas from a PFJet.