CMS 3D CMS Logo

PileupJetIdAlgo.cc
Go to the documentation of this file.
2 
11 
12 #include "TMatrixDSym.h"
13 #include "TMatrixDSymEigen.h"
14 
15 #include <utility>
16 
17 // ------------------------------------------------------------------------------------------
19 
20 // ------------------------------------------------------------------------------------------
21 
23  : cutBased_(ps.getParameter<bool>("cutBased")),
24  etaBinnedWeights_(false),
25  runMvas_(runMvas),
26  nEtaBins_(0),
27  label_(ps.getParameter<std::string>("label")),
28  mvacut_{},
29  rmsCut_{},
30  betaStarCut_{} {
31  std::vector<edm::FileInPath> tmvaEtaWeights;
32  std::vector<std::string> tmvaSpectators;
33  int version;
34 
35  if (!cutBased_) {
36  etaBinnedWeights_ = ps.getParameter<bool>("etaBinnedWeights");
37  if (etaBinnedWeights_) {
38  const std::vector<edm::ParameterSet>& trainings = ps.getParameter<std::vector<edm::ParameterSet> >("trainings");
39  nEtaBins_ = ps.getParameter<int>("nEtaBins");
40  for (int v = 0; v < nEtaBins_; v++) {
41  tmvaEtaWeights.push_back(trainings.at(v).getParameter<edm::FileInPath>("tmvaWeights"));
42  jEtaMin_.push_back(trainings.at(v).getParameter<double>("jEtaMin"));
43  jEtaMax_.push_back(trainings.at(v).getParameter<double>("jEtaMax"));
44  }
45  for (int v = 0; v < nEtaBins_; v++) {
46  tmvaEtaVariables_.push_back(trainings.at(v).getParameter<std::vector<std::string> >("tmvaVariables"));
47  }
48  } else {
49  tmvaVariables_ = ps.getParameter<std::vector<std::string> >("tmvaVariables");
50  }
51  tmvaMethod_ = ps.getParameter<std::string>("tmvaMethod");
52  tmvaSpectators = ps.getParameter<std::vector<std::string> >("tmvaSpectators");
53  version = ps.getParameter<int>("version");
54  } else {
55  version = USER;
56  }
57 
58  edm::ParameterSet jetConfig = ps.getParameter<edm::ParameterSet>("JetIdParams");
59  for (int i0 = 0; i0 < 3; i0++) {
60  std::string lCutType = "Tight";
62  lCutType = "Medium";
64  lCutType = "Loose";
65  int nCut = 1;
66  if (cutBased_)
67  nCut++;
68  for (int i1 = 0; i1 < nCut; i1++) {
69  std::string lFullCutType = lCutType;
70  if (cutBased_ && i1 == 0)
71  lFullCutType = "BetaStar" + lCutType;
72  if (cutBased_ && i1 == 1)
73  lFullCutType = "RMS" + lCutType;
74  std::vector<double> pt010 = jetConfig.getParameter<std::vector<double> >(("Pt010_" + lFullCutType).c_str());
75  std::vector<double> pt1020 = jetConfig.getParameter<std::vector<double> >(("Pt1020_" + lFullCutType).c_str());
76  std::vector<double> pt2030 = jetConfig.getParameter<std::vector<double> >(("Pt2030_" + lFullCutType).c_str());
77  std::vector<double> pt3040 = jetConfig.getParameter<std::vector<double> >(("Pt3040_" + lFullCutType).c_str());
78  std::vector<double> pt4050 = jetConfig.getParameter<std::vector<double> >(("Pt4050_" + lFullCutType).c_str());
79  if (!cutBased_) {
80  for (int i2 = 0; i2 < 4; i2++)
81  mvacut_[i0][0][i2] = pt010[i2];
82  for (int i2 = 0; i2 < 4; i2++)
83  mvacut_[i0][1][i2] = pt1020[i2];
84  for (int i2 = 0; i2 < 4; i2++)
85  mvacut_[i0][2][i2] = pt2030[i2];
86  for (int i2 = 0; i2 < 4; i2++)
87  mvacut_[i0][3][i2] = pt3040[i2];
88  for (int i2 = 0; i2 < 4; i2++)
89  mvacut_[i0][4][i2] = pt4050[i2];
90  }
91  if (cutBased_ && i1 == 0) {
92  for (int i2 = 0; i2 < 4; i2++)
93  betaStarCut_[i0][0][i2] = pt010[i2];
94  for (int i2 = 0; i2 < 4; i2++)
95  betaStarCut_[i0][1][i2] = pt1020[i2];
96  for (int i2 = 0; i2 < 4; i2++)
97  betaStarCut_[i0][2][i2] = pt2030[i2];
98  for (int i2 = 0; i2 < 4; i2++)
99  betaStarCut_[i0][3][i2] = pt3040[i2];
100  for (int i2 = 0; i2 < 4; i2++)
101  betaStarCut_[i0][4][i2] = pt4050[i2];
102  }
103  if (cutBased_ && i1 == 1) {
104  for (int i2 = 0; i2 < 4; i2++)
105  rmsCut_[i0][0][i2] = pt010[i2];
106  for (int i2 = 0; i2 < 4; i2++)
107  rmsCut_[i0][1][i2] = pt1020[i2];
108  for (int i2 = 0; i2 < 4; i2++)
109  rmsCut_[i0][2][i2] = pt2030[i2];
110  for (int i2 = 0; i2 < 4; i2++)
111  rmsCut_[i0][3][i2] = pt3040[i2];
112  for (int i2 = 0; i2 < 4; i2++)
113  rmsCut_[i0][4][i2] = pt4050[i2];
114  }
115  }
116  }
117 
118  if (!cutBased_) {
119  assert(tmvaMethod_.empty() || ((!tmvaVariables_.empty() || (!tmvaEtaVariables_.empty())) && version == USER));
120  }
121 
122  if ((!cutBased_) && (runMvas_)) {
123  if (etaBinnedWeights_) {
124  for (int v = 0; v < nEtaBins_; v++) {
125  etaReader_.push_back(createGBRForest(tmvaEtaWeights.at(v)));
126  }
127  } else {
128  reader_ = createGBRForest(ps.getParameter<std::string>("tmvaWeights"));
129  }
130  }
131 }
132 
134 
135 // ------------------------------------------------------------------------------------------
137 
138 // ------------------------------------------------------------------------------------------
139 void assign(const std::vector<float>& vec, float& a, float& b, float& c, float& d) {
140  size_t sz = vec.size();
141  a = (sz > 0 ? vec[0] : 0.);
142  b = (sz > 1 ? vec[1] : 0.);
143  c = (sz > 2 ? vec[2] : 0.);
144  d = (sz > 3 ? vec[3] : 0.);
145 }
146 // ------------------------------------------------------------------------------------------
147 void setPtEtaPhi(const reco::Candidate& p, float& pt, float& eta, float& phi) {
148  pt = p.pt();
149  eta = p.eta();
150  phi = p.phi();
151 }
152 
153 // ------------------------------------------------------------------------------------------
155 
156 // ------------------------------------------------------------------------------------------
157 
158 float PileupJetIdAlgo::getMVAval(const std::vector<std::string>& varList,
159  const std::unique_ptr<const GBRForest>& reader) {
160  std::vector<float> vars;
161  for (std::vector<std::string>::const_iterator it = varList.begin(); it != varList.end(); ++it) {
162  std::pair<float*, float> var = variables_.at(*it);
163  vars.push_back(*var.first);
164  }
165  return reader->GetClassifier(vars.data());
166 }
167 
169  if (cache_->cutBased()) {
171  internalId_.betaStarClassic_, internalId_.dR2Mean_, internalId_.nvtx_, internalId_.jetPt_, internalId_.jetEta_);
172  } else {
173  if (std::abs(internalId_.jetEta_) >= 5.0) {
174  internalId_.mva_ = -2.;
175  } else {
176  if (cache_->etaBinnedWeights()) {
177  if (std::abs(internalId_.jetEta_) > cache_->jEtaMax().at(cache_->nEtaBins() - 1)) {
178  internalId_.mva_ = -2.;
179  } else {
180  for (int v = 0; v < cache_->nEtaBins(); v++) {
181  if (std::abs(internalId_.jetEta_) >= cache_->jEtaMin().at(v) &&
182  std::abs(internalId_.jetEta_) < cache_->jEtaMax().at(v)) {
184  break;
185  }
186  }
187  }
188  } else {
190  }
191  }
193  }
194 }
195 
196 // ------------------------------------------------------------------------------------------
197 std::pair<int, int> PileupJetIdAlgo::getJetIdKey(float jetPt, float jetEta) {
198  int ptId = 0;
199  if (jetPt >= 10 && jetPt < 20)
200  ptId = 1;
201  if (jetPt >= 20 && jetPt < 30)
202  ptId = 2;
203  if (jetPt >= 30 && jetPt < 40)
204  ptId = 3;
205  if (jetPt >= 40)
206  ptId = 4;
207 
208  int etaId = 0;
209  if (std::abs(jetEta) >= 2.5 && std::abs(jetEta) < 2.75)
210  etaId = 1;
211  if (std::abs(jetEta) >= 2.75 && std::abs(jetEta) < 3.0)
212  etaId = 2;
213  if (std::abs(jetEta) >= 3.0 && std::abs(jetEta) < 5.0)
214  etaId = 3;
215 
216  return std::pair<int, int>(ptId, etaId);
217 }
218 // ------------------------------------------------------------------------------------------
219 int PileupJetIdAlgo::computeCutIDflag(float betaStarClassic, float dR2Mean, float nvtx, float jetPt, float jetEta) {
220  std::pair<int, int> jetIdKey = getJetIdKey(jetPt, jetEta);
221  float betaStarModified = betaStarClassic / log(nvtx - 0.64);
222  int idFlag(0);
223  if (betaStarModified < cache_->betaStarCut()[PileupJetIdentifier::kTight][jetIdKey.first][jetIdKey.second] &&
224  dR2Mean < cache_->rmsCut()[PileupJetIdentifier::kTight][jetIdKey.first][jetIdKey.second])
225  idFlag += 1 << PileupJetIdentifier::kTight;
226 
227  if (betaStarModified < cache_->betaStarCut()[PileupJetIdentifier::kMedium][jetIdKey.first][jetIdKey.second] &&
228  dR2Mean < cache_->rmsCut()[PileupJetIdentifier::kMedium][jetIdKey.first][jetIdKey.second])
229  idFlag += 1 << PileupJetIdentifier::kMedium;
230 
231  if (betaStarModified < cache_->betaStarCut()[PileupJetIdentifier::kLoose][jetIdKey.first][jetIdKey.second] &&
232  dR2Mean < cache_->rmsCut()[PileupJetIdentifier::kLoose][jetIdKey.first][jetIdKey.second])
233  idFlag += 1 << PileupJetIdentifier::kLoose;
234  return idFlag;
235 }
236 // ------------------------------------------------------------------------------------------
237 int PileupJetIdAlgo::computeIDflag(float mva, float jetPt, float jetEta) {
238  std::pair<int, int> jetIdKey = getJetIdKey(jetPt, jetEta);
239  return computeIDflag(mva, jetIdKey.first, jetIdKey.second);
240 }
241 
242 // ------------------------------------------------------------------------------------------
243 int PileupJetIdAlgo::computeIDflag(float mva, int ptId, int etaId) {
244  int idFlag(0);
245  if (mva > cache_->mvacut()[PileupJetIdentifier::kTight][ptId][etaId])
246  idFlag += 1 << PileupJetIdentifier::kTight;
247  if (mva > cache_->mvacut()[PileupJetIdentifier::kMedium][ptId][etaId])
248  idFlag += 1 << PileupJetIdentifier::kMedium;
249  if (mva > cache_->mvacut()[PileupJetIdentifier::kLoose][ptId][etaId])
250  idFlag += 1 << PileupJetIdentifier::kLoose;
251  return idFlag;
252 }
253 
254 // ------------------------------------------------------------------------------------------
256  runMva();
258 }
259 
260 // ------------------------------------------------------------------------------------------
262  float jec,
263  const reco::Vertex* vtx,
264  const reco::VertexCollection& allvtx,
265  double rho,
266  edm::ValueMap<float>& constituentWeights,
267  bool applyConstituentWeight) {
268  // initialize all variables to 0
269  resetVariables();
270 
271  // loop over constituents, accumulate sums and find leading candidates
272  const pat::Jet* patjet = dynamic_cast<const pat::Jet*>(jet);
273  const reco::PFJet* pfjet = dynamic_cast<const reco::PFJet*>(jet);
274  assert(patjet != nullptr || pfjet != nullptr);
275  if (patjet != nullptr && jec == 0.) { // if this is a pat jet and no jec has been passed take the jec from the object
276  jec = patjet->pt() / patjet->correctedJet(0).pt();
277  }
278  if (jec <= 0.) {
279  jec = 1.;
280  }
281 
282  const reco::Candidate *lLead = nullptr, *lSecond = nullptr, *lLeadNeut = nullptr, *lLeadEm = nullptr,
283  *lLeadCh = nullptr, *lTrail = nullptr;
284 
285  std::vector<float> frac, fracCh, fracEm, fracNeut;
286  float cones[] = {0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7};
287  size_t ncones = sizeof(cones) / sizeof(float);
288  float* coneFracs[] = {&internalId_.frac01_,
289  &internalId_.frac02_,
290  &internalId_.frac03_,
291  &internalId_.frac04_,
292  &internalId_.frac05_,
293  &internalId_.frac06_,
294  &internalId_.frac07_};
295  float* coneEmFracs[] = {&internalId_.emFrac01_,
302  float* coneNeutFracs[] = {&internalId_.neutFrac01_,
309  float* coneChFracs[] = {&internalId_.chFrac01_,
316  TMatrixDSym covMatrix(2);
317  covMatrix = 0.;
318  float jetPt = jet->pt() / jec; // use uncorrected pt for shape variables
319  float sumPt = 0., sumPt2 = 0., sumTkPt = 0., sumPtCh = 0, sumPtNe = 0;
320  float multNeut = 0.0;
321  float sumW2(0.0);
322  float sum_deta(0.0), sum_dphi(0.0);
323  float ave_deta(0.0), ave_dphi(0.0);
324  setPtEtaPhi(
325  *jet, internalId_.jetPt_, internalId_.jetEta_, internalId_.jetPhi_); // use corrected pt for jet kinematics
326  internalId_.jetM_ = jet->mass();
327  internalId_.nvtx_ = allvtx.size();
328  internalId_.rho_ = rho;
329 
330  float dRmin(1000);
331  float LeadcandWeight = 1.0;
332  float SecondcandWeight = 1.0;
333  float LeadNeutcandWeight = 1.0;
334  float LeadEmcandWeight = 1.0;
335  float LeadChcandWeight = 1.0;
336  float TrailcandWeight = 1.0;
337 
338  for (unsigned i = 0; i < jet->numberOfSourceCandidatePtrs(); ++i) {
339  reco::CandidatePtr pfJetConstituent = jet->sourceCandidatePtr(i);
340  const reco::Candidate* icand = pfJetConstituent.get();
341  const pat::PackedCandidate* lPack = dynamic_cast<const pat::PackedCandidate*>(icand);
342  const reco::PFCandidate* lPF = dynamic_cast<const reco::PFCandidate*>(icand);
343  bool isPacked = true;
344  if (lPack == nullptr) {
345  isPacked = false;
346  }
347 
348  float candWeight = 1.0;
349  if (applyConstituentWeight) { // PUPPI Jet weight should be pulled up from valuemap, not packed candidate
350  candWeight = constituentWeights[jet->sourceCandidatePtr(i)];
351  }
352  float candPt = (icand->pt()) * candWeight;
353  float candPtFrac = candPt / jetPt;
354  float candDr = reco::deltaR(*icand, *jet);
355  float candDeta = icand->eta() - jet->eta();
356  float candDphi = reco::deltaPhi(*icand, *jet);
357  float candPtDr = candPt * candDr;
358  size_t icone = std::lower_bound(&cones[0], &cones[ncones], candDr) - &cones[0];
359 
360  if (candDr < dRmin)
361  dRmin = candDr;
362 
363  // // all particles; PUPPI weights multiplied to leading and subleading constituent if it is for PUPPI
364  if (lLead == nullptr || candPt > (lLead->pt()) * LeadcandWeight) {
365  lSecond = lLead;
366  SecondcandWeight = LeadcandWeight;
367  lLead = icand;
368  if (applyConstituentWeight) {
369  LeadcandWeight = constituentWeights[jet->sourceCandidatePtr(i)];
370  }
371  } else if ((lSecond == nullptr || candPt > (lSecond->pt()) * SecondcandWeight) &&
372  (candPt < (lLead->pt()) * LeadcandWeight)) {
373  lSecond = icand;
374  if (applyConstituentWeight) {
375  SecondcandWeight = constituentWeights[jet->sourceCandidatePtr(i)];
376  }
377  }
378 
379  // // average shapes
380  internalId_.dRMean_ += candPtDr;
381  internalId_.dR2Mean_ += candPtDr * candPtDr;
382  covMatrix(0, 0) += candPt * candPt * candDeta * candDeta;
383  covMatrix(0, 1) += candPt * candPt * candDeta * candDphi;
384  covMatrix(1, 1) += candPt * candPt * candDphi * candDphi;
385  internalId_.ptD_ += candPt * candPt;
386  sumPt += candPt;
387  sumPt2 += candPt * candPt;
388 
389  // single most energetic candiates and jet shape profiles
390  frac.push_back(candPtFrac);
391 
392  if (icone < ncones) {
393  *coneFracs[icone] += candPt;
394  }
395 
396  // neutrals Neutral hadrons
397  if (abs(icand->pdgId()) == 130) {
398  if (lLeadNeut == nullptr || candPt > (lLeadNeut->pt()) * LeadNeutcandWeight) {
399  lLeadNeut = icand;
400  if (applyConstituentWeight) {
401  LeadNeutcandWeight = constituentWeights[jet->sourceCandidatePtr(i)];
402  }
403  }
404 
405  internalId_.dRMeanNeut_ += candPtDr;
406  fracNeut.push_back(candPtFrac);
407  if (icone < ncones) {
408  *coneNeutFracs[icone] += candPt;
409  }
410  internalId_.ptDNe_ += candPt * candPt;
411  sumPtNe += candPt;
412  multNeut += candWeight;
413  }
414 
415  // EM candidated photon
416  if (icand->pdgId() == 22) {
417  if (lLeadEm == nullptr || candPt > (lLeadEm->pt()) * LeadEmcandWeight) {
418  lLeadEm = icand;
419  if (applyConstituentWeight) {
420  LeadEmcandWeight = constituentWeights[jet->sourceCandidatePtr(i)];
421  }
422  }
423  internalId_.dRMeanEm_ += candPtDr;
424  fracEm.push_back(candPtFrac);
425  if (icone < ncones) {
426  *coneEmFracs[icone] += candPt;
427  }
428  internalId_.ptDNe_ += candPt * candPt;
429  sumPtNe += candPt;
430  multNeut += candWeight;
431  }
432  // hadrons and EM in HF
433  if ((abs(icand->pdgId()) == 1) || (abs(icand->pdgId()) == 2))
434  multNeut += candWeight;
435 
436  // Charged particles
437  if (icand->charge() != 0) {
438  if (lLeadCh == nullptr || candPt > (lLeadCh->pt()) * LeadChcandWeight) {
439  lLeadCh = icand;
440  if (applyConstituentWeight) {
441  LeadChcandWeight = constituentWeights[jet->sourceCandidatePtr(i)];
442  }
443  const reco::Track* pfTrk = icand->bestTrack();
444  if (lPF && std::abs(icand->pdgId()) == 13 && pfTrk == nullptr) {
445  reco::MuonRef lmuRef = lPF->muonRef();
446  if (lmuRef.isNonnull()) {
447  const reco::Muon& lmu = *lmuRef.get();
448  pfTrk = lmu.bestTrack();
449  edm::LogWarning("BadMuon")
450  << "Found a PFCandidate muon without a trackRef: falling back to Muon::bestTrack ";
451  }
452  }
453  if (pfTrk == nullptr) { //protection against empty pointers for the miniAOD case
454  //To handle the electron case
455  if (isPacked) {
456  internalId_.d0_ = std::abs(lPack->dxy(vtx->position()));
457  internalId_.dZ_ = std::abs(lPack->dz(vtx->position()));
458  } else if (lPF != nullptr) {
459  pfTrk = (lPF->trackRef().get() == nullptr) ? lPF->gsfTrackRef().get() : lPF->trackRef().get();
460  internalId_.d0_ = std::abs(pfTrk->dxy(vtx->position()));
461  internalId_.dZ_ = std::abs(pfTrk->dz(vtx->position()));
462  }
463  } else {
464  internalId_.d0_ = std::abs(pfTrk->dxy(vtx->position()));
465  internalId_.dZ_ = std::abs(pfTrk->dz(vtx->position()));
466  }
467  }
468  internalId_.dRMeanCh_ += candPtDr;
469  internalId_.ptDCh_ += candPt * candPt;
470  fracCh.push_back(candPtFrac);
471  if (icone < ncones) {
472  *coneChFracs[icone] += candPt;
473  }
474  sumPtCh += candPt;
475  }
476  // // beta and betastar
477  if (icand->charge() != 0) {
478  if (!isPacked) {
479  if (lPF->trackRef().isNonnull()) {
480  float tkpt = candPt;
481  sumTkPt += tkpt;
482  // 'classic' beta definition based on track-vertex association
483  bool inVtx0 = vtx->trackWeight(lPF->trackRef()) > 0;
484 
485  bool inAnyOther = false;
486  // alternative beta definition based on track-vertex distance of closest approach
487  double dZ0 = std::abs(lPF->trackRef()->dz(vtx->position()));
488  double dZ = dZ0;
489  for (reco::VertexCollection::const_iterator vi = allvtx.begin(); vi != allvtx.end(); ++vi) {
490  const reco::Vertex& iv = *vi;
491  if (iv.isFake() || iv.ndof() < 4) {
492  continue;
493  }
494  // the primary vertex may have been copied by the user: check identity by position
495  bool isVtx0 = (iv.position() - vtx->position()).r() < 0.02;
496  // 'classic' beta definition: check if the track is associated with
497  // any vertex other than the primary one
498  if (!isVtx0 && !inAnyOther) {
499  inAnyOther = vtx->trackWeight(lPF->trackRef()) <= 0;
500  }
501  // alternative beta: find closest vertex to the track
502  dZ = std::min(dZ, std::abs(lPF->trackRef()->dz(iv.position())));
503  }
504  // classic beta/betaStar
505  if (inVtx0 && !inAnyOther) {
506  internalId_.betaClassic_ += tkpt;
507  } else if (!inVtx0 && inAnyOther) {
508  internalId_.betaStarClassic_ += tkpt;
509  }
510  // alternative beta/betaStar
511  if (dZ0 < 0.2) {
512  internalId_.beta_ += tkpt;
513  } else if (dZ < 0.2) {
514  internalId_.betaStar_ += tkpt;
515  }
516  }
517  } else {
518  float tkpt = candPt;
519  sumTkPt += tkpt;
520  bool inVtx0 = false;
521  bool inVtxOther = false;
522  double dZ0 = 9999.;
523  double dZ_tmp = 9999.;
524  for (unsigned vtx_i = 0; vtx_i < allvtx.size(); vtx_i++) {
525  auto iv = allvtx[vtx_i];
526 
527  if (iv.isFake())
528  continue;
529 
530  // Match to vertex in case of copy as above
531  bool isVtx0 = (iv.position() - vtx->position()).r() < 0.02;
532 
533  if (isVtx0) {
534  if (lPack->fromPV(vtx_i) == pat::PackedCandidate::PVUsedInFit)
535  inVtx0 = true;
536  if (lPack->fromPV(vtx_i) == 0)
537  inVtxOther = true;
538  dZ0 = lPack->dz(iv.position());
539  }
540 
541  if (fabs(lPack->dz(iv.position())) < fabs(dZ_tmp)) {
542  dZ_tmp = lPack->dz(iv.position());
543  }
544  }
545  if (inVtx0) {
546  internalId_.betaClassic_ += tkpt;
547  } else if (inVtxOther) {
548  internalId_.betaStarClassic_ += tkpt;
549  }
550  if (fabs(dZ0) < 0.2) {
551  internalId_.beta_ += tkpt;
552  } else if (fabs(dZ_tmp) < 0.2) {
553  internalId_.betaStar_ += tkpt;
554  }
555  }
556  }
557 
558  // trailing candidate
559  if (lTrail == nullptr || candPt < (lTrail->pt()) * TrailcandWeight) {
560  lTrail = icand;
561  if (applyConstituentWeight) {
562  TrailcandWeight = constituentWeights[jet->sourceCandidatePtr(i)];
563  }
564  }
565 
566  // average for pull variable
567 
568  float weight2 = candPt * candPt;
569  sumW2 += weight2;
570  float deta = icand->eta() - jet->eta();
571  float dphi = reco::deltaPhi(*icand, *jet);
572  sum_deta += deta * weight2;
573  sum_dphi += dphi * weight2;
574  if (sumW2 > 0) {
575  ave_deta = sum_deta / sumW2;
576  ave_dphi = sum_dphi / sumW2;
577  }
578  }
579 
580  // // Finalize all variables
581  // Most of Below values are not used for puID variable generation at the moment, except lLeadCh Pt for JetRchg, so I assign that zero if there is no charged constituent.
582 
583  assert(!(lLead == nullptr));
584  internalId_.leadPt_ = lLead->pt() * LeadcandWeight;
585  internalId_.leadEta_ = lLead->eta();
586  internalId_.leadPhi_ = lLead->phi();
587 
588  if (lSecond != nullptr) {
589  internalId_.secondPt_ = lSecond->pt() * SecondcandWeight;
590  internalId_.secondEta_ = lSecond->eta();
591  internalId_.secondPhi_ = lSecond->phi();
592  } else {
593  internalId_.secondPt_ = 0.0;
596  }
597 
598  if (lLeadNeut != nullptr) {
599  internalId_.leadNeutPt_ = lLeadNeut->pt() * LeadNeutcandWeight;
600  internalId_.leadNeutEta_ = lLeadNeut->eta();
601  internalId_.leadNeutPhi_ = lLeadNeut->phi();
602  } else {
603  internalId_.leadNeutPt_ = 0.0;
606  }
607 
608  if (lLeadEm != nullptr) {
609  internalId_.leadEmPt_ = lLeadEm->pt() * LeadEmcandWeight;
610  internalId_.leadEmEta_ = lLeadEm->eta();
611  internalId_.leadEmPhi_ = lLeadEm->phi();
612  } else {
613  internalId_.leadEmPt_ = 0.0;
616  }
617 
618  if (lLeadCh != nullptr) {
619  internalId_.leadChPt_ = lLeadCh->pt() * LeadChcandWeight;
620  internalId_.leadChEta_ = lLeadCh->eta();
621  internalId_.leadChPhi_ = lLeadCh->phi();
622  } else {
623  internalId_.leadChPt_ = 0.0;
626  }
627 
628  if (patjet != nullptr) { // to enable running on MiniAOD slimmedJets
629  internalId_.nCharged_ = patjet->chargedMultiplicity();
630  internalId_.nNeutrals_ = patjet->neutralMultiplicity();
631  internalId_.chgEMfrac_ = patjet->chargedEmEnergy() / jet->energy();
632  internalId_.neuEMfrac_ = patjet->neutralEmEnergy() / jet->energy();
633  internalId_.chgHadrfrac_ = patjet->chargedHadronEnergy() / jet->energy();
634  internalId_.neuHadrfrac_ = patjet->neutralHadronEnergy() / jet->energy();
635  if (applyConstituentWeight)
636  internalId_.nNeutrals_ = multNeut;
637  } else {
638  internalId_.nCharged_ = pfjet->chargedMultiplicity();
639  internalId_.nNeutrals_ = pfjet->neutralMultiplicity();
640  internalId_.chgEMfrac_ = pfjet->chargedEmEnergy() / jet->energy();
641  internalId_.neuEMfrac_ = pfjet->neutralEmEnergy() / jet->energy();
642  internalId_.chgHadrfrac_ = pfjet->chargedHadronEnergy() / jet->energy();
643  internalId_.neuHadrfrac_ = pfjet->neutralHadronEnergy() / jet->energy();
644  if (applyConstituentWeight)
645  internalId_.nNeutrals_ = multNeut;
646  }
647  internalId_.nParticles_ = jet->nConstituents();
648 
650  float ddetaR_sum(0.0), ddphiR_sum(0.0), pull_tmp(0.0);
651  for (unsigned k = 0; k < jet->numberOfSourceCandidatePtrs(); k++) {
652  reco::CandidatePtr temp_pfJetConstituent = jet->sourceCandidatePtr(k);
653  // reco::CandidatePtr temp_weightpfJetConstituent = jet->sourceCandidatePtr(k);
654  const reco::Candidate* part = temp_pfJetConstituent.get();
655 
656  float candWeight = 1.0;
657 
658  if (applyConstituentWeight)
659  candWeight = constituentWeights[jet->sourceCandidatePtr(k)];
660 
661  float weight = candWeight * (part->pt()) * candWeight * (part->pt());
662  float deta = part->eta() - jet->eta();
663  float dphi = reco::deltaPhi(*part, *jet);
664  float ddeta, ddphi, ddR;
665  ddeta = deta - ave_deta;
666  ddphi = dphi - ave_dphi;
667  ddR = sqrt(ddeta * ddeta + ddphi * ddphi);
668  ddetaR_sum += ddR * ddeta * weight;
669  ddphiR_sum += ddR * ddphi * weight;
670  }
671  if (sumW2 > 0) {
672  float ddetaR_ave = ddetaR_sum / sumW2;
673  float ddphiR_ave = ddphiR_sum / sumW2;
674  pull_tmp = sqrt(ddetaR_ave * ddetaR_ave + ddphiR_ave * ddphiR_ave);
675  }
676  internalId_.pull_ = pull_tmp;
678 
679  std::sort(frac.begin(), frac.end(), std::greater<float>());
680  std::sort(fracCh.begin(), fracCh.end(), std::greater<float>());
681  std::sort(fracEm.begin(), fracEm.end(), std::greater<float>());
682  std::sort(fracNeut.begin(), fracNeut.end(), std::greater<float>());
684  assign(
686  assign(
688  assign(fracNeut,
693 
694  covMatrix(0, 0) /= sumPt2;
695  covMatrix(0, 1) /= sumPt2;
696  covMatrix(1, 1) /= sumPt2;
697  covMatrix(1, 0) = covMatrix(0, 1);
698  internalId_.etaW_ = sqrt(covMatrix(0, 0));
699  internalId_.phiW_ = sqrt(covMatrix(1, 1));
701  TVectorD eigVals(2);
702  eigVals = TMatrixDSymEigen(covMatrix).GetEigenValues();
703  internalId_.majW_ = sqrt(std::abs(eigVals(0)));
704  internalId_.minW_ = sqrt(std::abs(eigVals(1)));
705  if (internalId_.majW_ < internalId_.minW_) {
706  std::swap(internalId_.majW_, internalId_.minW_);
707  }
708 
710  if (lSecond != nullptr) {
711  internalId_.dRLead2nd_ = reco::deltaR(*jet, *lSecond);
712  }
713  internalId_.dRMean_ /= jetPt;
717  internalId_.dR2Mean_ /= sumPt2;
718 
719  for (size_t ic = 0; ic < ncones; ++ic) {
720  *coneFracs[ic] /= jetPt;
721  *coneEmFracs[ic] /= jetPt;
722  *coneNeutFracs[ic] /= jetPt;
723  *coneChFracs[ic] /= jetPt;
724  }
725  //http://jets.physics.harvard.edu/qvg/
726  double ptMean = sumPt / internalId_.nParticles_;
727  double ptRMS = 0;
728  for (unsigned int i0 = 0; i0 < frac.size(); i0++) {
729  ptRMS += (frac[i0] - ptMean) * (frac[i0] - ptMean);
730  }
731  ptRMS /= internalId_.nParticles_;
732  ptRMS = sqrt(ptRMS);
733 
734  internalId_.ptMean_ = ptMean;
735  internalId_.ptRMS_ = ptRMS / jetPt;
736  internalId_.pt2A_ = sqrt(internalId_.ptD_ / internalId_.nParticles_) / jetPt;
737  internalId_.ptD_ = sqrt(internalId_.ptD_) / sumPt;
741  internalId_.sumChPt_ = sumPtCh;
742  internalId_.sumNePt_ = sumPtNe;
743 
744  internalId_.jetR_ = (lLead->pt()) * LeadcandWeight / sumPt;
745  if (lLeadCh != nullptr) {
746  internalId_.jetRchg_ = (lLeadCh->pt()) * LeadChcandWeight / sumPt;
747  } else {
748  internalId_.jetRchg_ = 0;
749  }
750 
751  internalId_.dRMatch_ = dRmin;
752 
753  if (sumTkPt != 0.) {
754  internalId_.beta_ /= sumTkPt;
755  internalId_.betaStar_ /= sumTkPt;
756  internalId_.betaClassic_ /= sumTkPt;
757  internalId_.betaStarClassic_ /= sumTkPt;
758  } else {
759  assert(internalId_.beta_ == 0. && internalId_.betaStar_ == 0. && internalId_.betaClassic_ == 0. &&
760  internalId_.betaStarClassic_ == 0.);
761  }
762 
763  if (cache_->runMvas()) {
764  runMva();
765  }
766 
768 }
769 
770 // ------------------------------------------------------------------------------------------
772  std::stringstream out;
773  for (variables_list_t::const_iterator it = variables_.begin(); it != variables_.end(); ++it) {
774  out << std::setw(15) << it->first << std::setw(3) << "=" << std::setw(5) << *it->second.first << " ("
775  << std::setw(5) << it->second.second << ")" << std::endl;
776  }
777  return out.str();
778 }
779 
780 // ------------------------------------------------------------------------------------------
782  internalId_.idFlag_ = 0;
783  for (variables_list_t::iterator it = variables_.begin(); it != variables_.end(); ++it) {
784  *it->second.first = it->second.second;
785  }
786 }
787 
788 // ------------------------------------------------------------------------------------------
789 #define INIT_VARIABLE(NAME, TMVANAME, VAL) \
790  internalId_.NAME##_ = VAL; \
791  variables_[#NAME] = std::make_pair(&internalId_.NAME##_, VAL);
792 
793 // ------------------------------------------------------------------------------------------
795  internalId_.idFlag_ = 0;
796  INIT_VARIABLE(mva, "", -100.);
797  //INIT_VARIABLE(jetPt , "jspt_1", 0.);
798  //INIT_VARIABLE(jetEta , "jseta_1", large_val);
799  INIT_VARIABLE(jetPt, "", 0.);
801  INIT_VARIABLE(jetPhi, "jsphi_1", large_val);
802  INIT_VARIABLE(jetM, "jm_1", 0.);
803 
804  INIT_VARIABLE(nCharged, "", 0.);
805  INIT_VARIABLE(nNeutrals, "", 0.);
806 
807  INIT_VARIABLE(chgEMfrac, "", 0.);
808  INIT_VARIABLE(neuEMfrac, "", 0.);
809  INIT_VARIABLE(chgHadrfrac, "", 0.);
810  INIT_VARIABLE(neuHadrfrac, "", 0.);
811 
812  INIT_VARIABLE(d0, "jd0_1", -1000.);
813  INIT_VARIABLE(dZ, "jdZ_1", -1000.);
814  //INIT_VARIABLE(nParticles , "npart_1" , 0.);
815  INIT_VARIABLE(nParticles, "", 0.);
816 
817  INIT_VARIABLE(leadPt, "lpt_1", 0.);
818  INIT_VARIABLE(leadEta, "leta_1", large_val);
819  INIT_VARIABLE(leadPhi, "lphi_1", large_val);
820  INIT_VARIABLE(secondPt, "spt_1", 0.);
821  INIT_VARIABLE(secondEta, "seta_1", large_val);
822  INIT_VARIABLE(secondPhi, "sphi_1", large_val);
823  INIT_VARIABLE(leadNeutPt, "lnept_1", 0.);
824  INIT_VARIABLE(leadNeutEta, "lneeta_1", large_val);
825  INIT_VARIABLE(leadNeutPhi, "lnephi_1", large_val);
826  INIT_VARIABLE(leadEmPt, "lempt_1", 0.);
827  INIT_VARIABLE(leadEmEta, "lemeta_1", large_val);
828  INIT_VARIABLE(leadEmPhi, "lemphi_1", large_val);
829  INIT_VARIABLE(leadChPt, "lchpt_1", 0.);
830  INIT_VARIABLE(leadChEta, "lcheta_1", large_val);
831  INIT_VARIABLE(leadChPhi, "lchphi_1", large_val);
832  INIT_VARIABLE(leadFrac, "lLfr_1", 0.);
833 
834  INIT_VARIABLE(dRLeadCent, "drlc_1", 0.);
835  INIT_VARIABLE(dRLead2nd, "drls_1", 0.);
836  INIT_VARIABLE(dRMean, "drm_1", 0.);
837  INIT_VARIABLE(dRMean, "", 0.);
838  INIT_VARIABLE(pull, "", 0.);
839  INIT_VARIABLE(dRMeanNeut, "drmne_1", 0.);
840  INIT_VARIABLE(dRMeanEm, "drem_1", 0.);
841  INIT_VARIABLE(dRMeanCh, "drch_1", 0.);
842  INIT_VARIABLE(dR2Mean, "", 0.);
843 
844  INIT_VARIABLE(ptD, "", 0.);
845  INIT_VARIABLE(ptMean, "", 0.);
846  INIT_VARIABLE(ptRMS, "", 0.);
847  INIT_VARIABLE(pt2A, "", 0.);
848  INIT_VARIABLE(ptDCh, "", 0.);
849  INIT_VARIABLE(ptDNe, "", 0.);
850  INIT_VARIABLE(sumPt, "", 0.);
851  INIT_VARIABLE(sumChPt, "", 0.);
852  INIT_VARIABLE(sumNePt, "", 0.);
853 
854  INIT_VARIABLE(secondFrac, "", 0.);
855  INIT_VARIABLE(thirdFrac, "", 0.);
856  INIT_VARIABLE(fourthFrac, "", 0.);
857 
858  INIT_VARIABLE(leadChFrac, "", 0.);
859  INIT_VARIABLE(secondChFrac, "", 0.);
860  INIT_VARIABLE(thirdChFrac, "", 0.);
861  INIT_VARIABLE(fourthChFrac, "", 0.);
862 
863  INIT_VARIABLE(leadNeutFrac, "", 0.);
864  INIT_VARIABLE(secondNeutFrac, "", 0.);
865  INIT_VARIABLE(thirdNeutFrac, "", 0.);
866  INIT_VARIABLE(fourthNeutFrac, "", 0.);
867 
868  INIT_VARIABLE(leadEmFrac, "", 0.);
869  INIT_VARIABLE(secondEmFrac, "", 0.);
870  INIT_VARIABLE(thirdEmFrac, "", 0.);
871  INIT_VARIABLE(fourthEmFrac, "", 0.);
872 
873  INIT_VARIABLE(jetW, "", 1.);
874  INIT_VARIABLE(etaW, "", 1.);
875  INIT_VARIABLE(phiW, "", 1.);
876 
877  INIT_VARIABLE(majW, "", 1.);
878  INIT_VARIABLE(minW, "", 1.);
879 
880  INIT_VARIABLE(frac01, "", 0.);
881  INIT_VARIABLE(frac02, "", 0.);
882  INIT_VARIABLE(frac03, "", 0.);
883  INIT_VARIABLE(frac04, "", 0.);
884  INIT_VARIABLE(frac05, "", 0.);
885  INIT_VARIABLE(frac06, "", 0.);
886  INIT_VARIABLE(frac07, "", 0.);
887 
888  INIT_VARIABLE(chFrac01, "", 0.);
889  INIT_VARIABLE(chFrac02, "", 0.);
890  INIT_VARIABLE(chFrac03, "", 0.);
891  INIT_VARIABLE(chFrac04, "", 0.);
892  INIT_VARIABLE(chFrac05, "", 0.);
893  INIT_VARIABLE(chFrac06, "", 0.);
894  INIT_VARIABLE(chFrac07, "", 0.);
895 
896  INIT_VARIABLE(neutFrac01, "", 0.);
897  INIT_VARIABLE(neutFrac02, "", 0.);
898  INIT_VARIABLE(neutFrac03, "", 0.);
899  INIT_VARIABLE(neutFrac04, "", 0.);
900  INIT_VARIABLE(neutFrac05, "", 0.);
901  INIT_VARIABLE(neutFrac06, "", 0.);
902  INIT_VARIABLE(neutFrac07, "", 0.);
903 
904  INIT_VARIABLE(emFrac01, "", 0.);
905  INIT_VARIABLE(emFrac02, "", 0.);
906  INIT_VARIABLE(emFrac03, "", 0.);
907  INIT_VARIABLE(emFrac04, "", 0.);
908  INIT_VARIABLE(emFrac05, "", 0.);
909  INIT_VARIABLE(emFrac06, "", 0.);
910  INIT_VARIABLE(emFrac07, "", 0.);
911 
912  INIT_VARIABLE(beta, "", 0.);
913  INIT_VARIABLE(betaStar, "", 0.);
914  INIT_VARIABLE(betaClassic, "", 0.);
915  INIT_VARIABLE(betaStarClassic, "", 0.);
916 
917  INIT_VARIABLE(nvtx, "", 0.);
918  INIT_VARIABLE(rho, "", 0.);
919  INIT_VARIABLE(nTrueInt, "", 0.);
920 
921  INIT_VARIABLE(jetR, "", 0.);
922  INIT_VARIABLE(jetRchg, "", 0.);
923  INIT_VARIABLE(dRMatch, "", 0.);
924 }
925 
926 #undef INIT_VARIABLE
const float large_val
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
void set(const PileupJetIdentifier &)
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:467
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< double > const & jEtaMin() const
AlgoGBRForestsAndConstants const * cache_
#define INIT_VARIABLE(NAME, TMVANAME, VAL)
AlgoGBRForestsAndConstants(edm::ParameterSet const &, bool runMvas)
int32_t *__restrict__ iv
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:99
double pt() const final
transverse momentum
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:152
int chargedMultiplicity() const
chargedMultiplicity
Definition: Jet.h:740
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:139
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:147
virtual double pt() const =0
transverse momentum
Base class for all types of Jets.
Definition: Jet.h:20
std::string dumpVariables() const
std::pair< int, int > getJetIdKey(float jetPt, float jetEta)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Definition: weight.py:1
const PVAssoc fromPV(size_t ipv=0) const
int neutralMultiplicity() const
neutralMultiplicity
Definition: Jet.h:464
reader
Definition: DQM.py:105
int computeCutIDflag(float betaStarClassic, float dR2Mean, float nvtx, float jetPt, float jetEta)
const Track * bestTrack() const override
Definition: Muon.h:58
int neutralMultiplicity() const
neutralMultiplicity
Definition: PFJet.h:154
assert(be >=bs)
Jets made from PFObjects.
Definition: PFJet.h:20
std::vector< std::vector< std::string > > const & tmvaEtaVariables() const
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: Jet.h:694
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:623
PileupJetIdentifier internalId_
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho, edm::ValueMap< float > &constituentWeights, bool applyConstituentWeight)
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::string > const & tmvaVariables() const
void assign(const std::vector< float > &vec, float &a, float &b, float &c, float &d)
virtual int charge() const =0
electric charge
Jet correctedJet(const std::string &level, const std::string &flavor="none", const std::string &set="") const
float chargedEmEnergy() const
chargedEmEnergy
Definition: Jet.h:712
reco::MuonRef muonRef() const
Definition: PFCandidate.cc:448
d
Definition: ztail.py:151
int computeIDflag(float mva, float jetPt, float jetEta)
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
virtual int pdgId() const =0
PDG identifier.
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: Jet.h:704
void setPtEtaPhi(const reco::Candidate &p, float &pt, float &eta, float &phi)
static constexpr float d0
std::vector< std::unique_ptr< const GBRForest > > const & etaReader() const
part
Definition: HCALResponse.h:20
PileupJetIdentifier computeMva()
double b
Definition: hdecay.h:118
Analysis-level calorimeter jet class.
Definition: Jet.h:77
std::unique_ptr< const GBRForest > const & reader() const
float neutralEmEnergy() const
neutralEmEnergy
Definition: Jet.h:722
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
double a
Definition: hdecay.h:119
def cache(function)
Definition: utilities.py:3
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:95
virtual const Track * bestTrack() const
Definition: Candidate.h:268
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
variables_list_t variables_
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:435
virtual float dxy() const
dxy with respect to the PV ref
vars
Definition: DeepTauId.cc:30
virtual float dz(size_t ipv=0) const
dz with respect to the PV[ipv]
Log< level::Warning, false > LogWarning
std::vector< double > const & jEtaMax() const
float getMVAval(const std::vector< std::string > &, const std::unique_ptr< const GBRForest > &)
PileupJetIdAlgo(AlgoGBRForestsAndConstants const *cache)
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:609