CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Static Private Attributes
MVAJetPuId Class Reference

#include <MVAJetPuId.h>

Public Types

typedef std::map< std::string, std::pair< float *, float > > variables_list_t
 
enum  version_t { USER = -1, CATEv0 = 0 }
 

Public Member Functions

int computeCutIDflag (float betaStarClassic, float dR2Mean, float nvtx, float jetPt, float jetEta)
 
int computeIDflag (float mva, float jetPt, float jetEta)
 
int computeIDflag (float mva, int ptId, int etaId)
 
PileupJetIdentifier computeIdVariables (const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho, bool calculateMva=false)
 
PileupJetIdentifier computeMva ()
 
std::string dumpVariables () const
 
std::pair< int, int > getJetIdKey (float jetPt, float jetEta)
 
const variables_list_tgetVariables () const
 
const std::string method () const
 
 MVAJetPuId (int version=CATEv0, const std::string &tmvaWeight="", const std::string &tmvaMethod="", Float_t impactParTkThreshod_=1., const std::vector< std::string > &tmvaVariables=std::vector< std::string >())
 
 MVAJetPuId (const edm::ParameterSet &ps)
 
void set (const PileupJetIdentifier &)
 
 ~MVAJetPuId ()
 

Protected Member Functions

void bookReader ()
 
void initVariables ()
 
void resetVariables ()
 
void runMva ()
 
void setup ()
 

Protected Attributes

Float_t betaStarCut_ [NWPs][NEtas][NPts]
 
bool cutBased_
 
Float_t impactParTkThreshod_
 
PileupJetIdentifier internalId_
 
Float_t mvacut_ [NWPs][NEtas][NPts]
 
TMVA::Reader * reader_
 
Float_t rmsCut_ [NWPs][NEtas][NPts]
 
std::string tmvaMethod_
 
std::map< std::string, std::string > tmvaNames_
 
std::vector< std::string > tmvaSpectators_
 
std::vector< std::string > tmvaVariables_
 
std::string tmvaWeights_
 
variables_list_t variables_
 
Int_t version_
 

Static Private Attributes

static constexpr int NEtas = 4
 
static constexpr int NPts = 4
 
static constexpr int NWPs = 3
 

Detailed Description

Definition at line 16 of file MVAJetPuId.h.

Member Typedef Documentation

◆ variables_list_t

typedef std::map<std::string, std::pair<float *, float> > MVAJetPuId::variables_list_t

Definition at line 45 of file MVAJetPuId.h.

Member Enumeration Documentation

◆ version_t

Enumerator
USER 
CATEv0 

Definition at line 22 of file MVAJetPuId.h.

Constructor & Destructor Documentation

◆ MVAJetPuId() [1/2]

MVAJetPuId::MVAJetPuId ( int  version = CATEv0,
const std::string &  tmvaWeight = "",
const std::string &  tmvaMethod = "",
Float_t  impactParTkThreshod_ = 1.,
const std::vector< std::string > &  tmvaVariables = std::vector<std::string>() 
)

Definition at line 50 of file MVAJetPuId.cc.

References impactParTkThreshod_, reader_, setup(), hltPUIdAlgo_cff::tmvaMethod, tmvaMethod_, hltPUIdAlgo_cff::tmvaVariables, tmvaVariables_, hltPUIdAlgo_cff::tmvaWeights, tmvaWeights_, BeamSplash_cfg::version, and version_.

54  {
55  impactParTkThreshod_ = impactParTkThreshod;
59  version_ = version;
60 
61  reader_ = nullptr;
62 
63  setup();
64 }
std::vector< std::string > tmvaVariables_
Definition: MVAJetPuId.h:66
Int_t version_
Definition: MVAJetPuId.h:70
TMVA::Reader * reader_
Definition: MVAJetPuId.h:64
std::string tmvaWeights_
Definition: MVAJetPuId.h:65
Float_t impactParTkThreshod_
Definition: MVAJetPuId.h:71
void setup()
Definition: MVAJetPuId.cc:66
std::string tmvaMethod_
Definition: MVAJetPuId.h:65

◆ MVAJetPuId() [2/2]

MVAJetPuId::MVAJetPuId ( const edm::ParameterSet ps)

Definition at line 16 of file MVAJetPuId.cc.

References edm::FileInPath::fullPath(), edm::ParameterSet::getParameter(), testProducerWithPsetDescEmpty_cfi::i1, testProducerWithPsetDescEmpty_cfi::i2, impactParTkThreshod_, PileupJetIdentifier::kLoose, PileupJetIdentifier::kMedium, mvacut_, NPts, NWPs, reader_, setup(), AlCaHLTBitMon_QueryRunRegistry::string, tmvaMethod_, tmvaSpectators_, tmvaVariables_, tmvaWeights_, and version_.

16  {
18 
20  tmvaMethod_ = ps.getParameter<std::string>("tmvaMethod");
21  tmvaVariables_ = ps.getParameter<std::vector<std::string> >("tmvaVariables");
22  tmvaSpectators_ = ps.getParameter<std::vector<std::string> >("tmvaSpectators");
23  version_ = ps.getParameter<int>("version");
24  reader_ = nullptr;
25  edm::ParameterSet jetConfig = ps.getParameter<edm::ParameterSet>("JetIdParams");
26  for (int i0 = 0; i0 < NWPs; i0++) {
27  std::string lCutType = "Tight";
29  lCutType = "Medium";
31  lCutType = "Loose";
32  for (int i1 = 0; i1 < 1; i1++) {
33  std::vector<double> pt010 = jetConfig.getParameter<std::vector<double> >(("Pt010_" + lCutType).c_str());
34  std::vector<double> pt1020 = jetConfig.getParameter<std::vector<double> >(("Pt1020_" + lCutType).c_str());
35  std::vector<double> pt2030 = jetConfig.getParameter<std::vector<double> >(("Pt2030_" + lCutType).c_str());
36  std::vector<double> pt3050 = jetConfig.getParameter<std::vector<double> >(("Pt3050_" + lCutType).c_str());
37  for (int i2 = 0; i2 < NPts; i2++)
38  mvacut_[i0][0][i2] = pt010[i2];
39  for (int i2 = 0; i2 < NPts; i2++)
40  mvacut_[i0][1][i2] = pt1020[i2];
41  for (int i2 = 0; i2 < NPts; i2++)
42  mvacut_[i0][2][i2] = pt2030[i2];
43  for (int i2 = 0; i2 < NPts; i2++)
44  mvacut_[i0][3][i2] = pt3050[i2];
45  }
46  }
47  setup();
48 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< std::string > tmvaVariables_
Definition: MVAJetPuId.h:66
std::string fullPath() const
Definition: FileInPath.cc:161
Float_t mvacut_[NWPs][NEtas][NPts]
Definition: MVAJetPuId.h:73
Int_t version_
Definition: MVAJetPuId.h:70
std::vector< std::string > tmvaSpectators_
Definition: MVAJetPuId.h:67
static constexpr int NPts
Definition: MVAJetPuId.h:18
TMVA::Reader * reader_
Definition: MVAJetPuId.h:64
std::string tmvaWeights_
Definition: MVAJetPuId.h:65
Float_t impactParTkThreshod_
Definition: MVAJetPuId.h:71
void setup()
Definition: MVAJetPuId.cc:66
static constexpr int NWPs
Definition: MVAJetPuId.h:17
std::string tmvaMethod_
Definition: MVAJetPuId.h:65

◆ ~MVAJetPuId()

MVAJetPuId::~MVAJetPuId ( )

Definition at line 105 of file MVAJetPuId.cc.

References reader_.

105  {
106  if (reader_) {
107  delete reader_;
108  }
109 }
TMVA::Reader * reader_
Definition: MVAJetPuId.h:64

Member Function Documentation

◆ bookReader()

void MVAJetPuId::bookReader ( )
protected

Definition at line 125 of file MVAJetPuId.cc.

References cms::cuda::assert(), relativeConstraints::empty, first, reco::details::loadTMVAWeights(), reader_, tmvaMethod_, tmvaNames_, tmvaSpectators_, tmvaVariables_, tmvaWeights_, and variables_.

Referenced by runMva().

125  {
126  reader_ = new TMVA::Reader("!Color:Silent");
127  assert(!tmvaMethod_.empty() && !tmvaWeights_.empty());
128  for (std::vector<std::string>::iterator it = tmvaVariables_.begin(); it != tmvaVariables_.end(); ++it) {
129  if (tmvaNames_[*it].empty()) {
130  tmvaNames_[*it] = *it;
131  }
132  reader_->AddVariable(*it, variables_[tmvaNames_[*it]].first);
133  }
134  for (std::vector<std::string>::iterator it = tmvaSpectators_.begin(); it != tmvaSpectators_.end(); ++it) {
135  if (tmvaNames_[*it].empty()) {
136  tmvaNames_[*it] = *it;
137  }
138  reader_->AddSpectator(*it, variables_[tmvaNames_[*it]].first);
139  }
141 }
std::vector< std::string > tmvaVariables_
Definition: MVAJetPuId.h:66
std::vector< std::string > tmvaSpectators_
Definition: MVAJetPuId.h:67
variables_list_t variables_
Definition: MVAJetPuId.h:62
assert(be >=bs)
TMVA::Reader * reader_
Definition: MVAJetPuId.h:64
TMVA::IMethod * loadTMVAWeights(TMVA::Reader *reader, const std::string &method, const std::string &weightFile, bool verbose=false)
std::string tmvaWeights_
Definition: MVAJetPuId.h:65
std::map< std::string, std::string > tmvaNames_
Definition: MVAJetPuId.h:68
std::string tmvaMethod_
Definition: MVAJetPuId.h:65

◆ computeCutIDflag()

int MVAJetPuId::computeCutIDflag ( float  betaStarClassic,
float  dR2Mean,
float  nvtx,
float  jetPt,
float  jetEta 
)

◆ computeIDflag() [1/2]

int MVAJetPuId::computeIDflag ( float  mva,
float  jetPt,
float  jetEta 
)

Definition at line 175 of file MVAJetPuId.cc.

References getJetIdKey(), reco::btau::jetEta, reco::btau::jetPt, and beam_dqm_sourceclient-live_cfg::mva.

Referenced by runMva().

175  {
176  std::pair<int, int> jetIdKey = getJetIdKey(jetPt, jetEta);
177  return computeIDflag(mva, jetIdKey.first, jetIdKey.second);
178 }
std::pair< int, int > getJetIdKey(float jetPt, float jetEta)
Definition: MVAJetPuId.cc:156
int computeIDflag(float mva, float jetPt, float jetEta)
Definition: MVAJetPuId.cc:175

◆ computeIDflag() [2/2]

int MVAJetPuId::computeIDflag ( float  mva,
int  ptId,
int  etaId 
)

◆ computeIdVariables()

PileupJetIdentifier MVAJetPuId::computeIdVariables ( const reco::Jet jet,
float  jec,
const reco::Vertex vtx,
const reco::VertexCollection allvtx,
double  rho,
bool  calculateMva = false 
)

Definition at line 196 of file MVAJetPuId.cc.

References a, cms::cuda::assert(), Assign(), b, c, reco::PFJet::chargedEmEnergy(), reco::PFJet::chargedHadronEnergy(), reco::PFJet::chargedMultiplicity(), PileupJetIdentifier::chgEMfrac_, PileupJetIdentifier::chgHadrfrac_, PileupJetIdentifier::d0_, dumpMFGeometry_cfg::delta, reco::deltaPhi(), reco::deltaR(), PileupJetIdentifier::dRLead2nd_, PileupJetIdentifier::dRMeanEm_, PileupJetIdentifier::dRMeanNeut_, spr::find(), PileupJetIdentifier::fourthFrac_, DivergingColor::frac, reco::PFCandidate::gamma, reco::PFJet::getPFConstituent(), reco::PFJet::getPFConstituents(), reco::PFCandidate::h0, mps_fire::i, internalId_, edm::Ptr< T >::isAvailable(), edm::Ref< C, T, F >::isAvailable(), edm::Ptr< T >::isNonnull(), edm::Ref< C, T, F >::isNonnull(), edm::Ptr< T >::isNull(), gpuVertexFinder::iv, jetMETDQMOfflineSource_cff::jec, metsig::jet, PileupJetIdentifier::jetM_, PileupJetIdentifier::jetPhi_, reco::btau::jetPt, PileupJetIdentifier::leadChEta_, PileupJetIdentifier::leadChPhi_, PileupJetIdentifier::leadChPt_, PileupJetIdentifier::leadEmEta_, PileupJetIdentifier::leadEmPhi_, PileupJetIdentifier::leadEmPt_, PileupJetIdentifier::leadEta_, PileupJetIdentifier::leadFrac_, PileupJetIdentifier::leadNeutEta_, PileupJetIdentifier::leadNeutPhi_, PileupJetIdentifier::leadNeutPt_, PileupJetIdentifier::leadPhi_, PileupJetIdentifier::leadPt_, pfDeepBoostedJetPreprocessParams_cfi::lower_bound, SiStripPI::min, PileupJetIdentifier::neuEMfrac_, PileupJetIdentifier::neuHadrfrac_, reco::PFJet::neutralEmEnergy(), reco::PFJet::neutralHadronEnergy(), reco::PFJet::neutralMultiplicity(), DiDispStaMuonMonitor_cfi::pt, PileupJetIdentifier::pt2A_, PileupJetIdentifier::ptMean_, PileupJetIdentifier::ptRMS_, alignCSCRings::r, resetVariables(), rho, runMva(), PileupJetIdentifier::secondEta_, PileupJetIdentifier::secondFrac_, PileupJetIdentifier::secondPhi_, PileupJetIdentifier::secondPt_, SetPtEtaPhi(), jetUpdater_cfi::sort, mathSSE::sqrt(), PileupJetIdentifier::sumChPt_, PileupJetIdentifier::sumNePt_, TtFullHadEvtBuilder_cfi::sumPt, PileupJetIdentifier::sumPt_, std::swap(), PileupJetIdentifier::thirdFrac_, extraflags_cff::vtx, and mps_merge::weight.

Referenced by MVAJetPuIdProducer::produce().

201  {
202  typedef std::vector<reco::PFCandidatePtr> constituents_type;
203  typedef std::vector<reco::PFCandidatePtr>::iterator constituents_iterator;
204 
205  resetVariables();
206 
207  const reco::PFJet *pfjet = dynamic_cast<const reco::PFJet *>(jet);
208 
209  if (jec < 0.) {
210  jec = 1.;
211  }
212 
213  constituents_type constituents = pfjet->getPFConstituents();
214 
215  reco::PFCandidatePtr lLead, lSecond, lLeadNeut, lLeadEm, lLeadCh, lTrail;
216  std::vector<float> frac, fracCh, fracEm, fracNeut;
217  constexpr int ncones = 4;
218  std::array<float, ncones> cones{{0.1, 0.2, 0.3, 0.4}};
219  std::array<float *, ncones> coneFracs{
220  {&internalId_.frac01_, &internalId_.frac02_, &internalId_.frac03_, &internalId_.frac04_}};
221  TMatrixDSym covMatrix(2);
222  covMatrix = 0.;
223 
224  reco::TrackRef impactTrack;
225  float jetPt = jet->pt() / jec; // use uncorrected pt for shape variables
226  float sumPt = 0., sumPt2 = 0., sumTkPt = 0., sumPtCh = 0, sumPtNe = 0;
227  float sum_deta = 0;
228  float sum_dphi = 0;
229  float sum_deta2 = 0;
230  float sum_detadphi = 0;
231  float sum_dphi2 = 0;
232  SetPtEtaPhi(
233  *jet, internalId_.jetPt_, internalId_.jetEta_, internalId_.jetPhi_); // use corrected pt for jet kinematics
234  internalId_.jetM_ = jet->mass();
235  internalId_.rho_ = rho; //allvtx.size();
236  for (constituents_iterator it = constituents.begin(); it != constituents.end(); ++it) {
237  reco::PFCandidatePtr &icand = *it;
238  float candPt = icand->pt();
239  float candPtFrac = candPt / jetPt;
240  float candDr = reco::deltaR(**it, *jet);
241  float candDeta = fabs((*it)->eta() - jet->eta());
242  float candDphi = reco::deltaPhi(**it, *jet);
243  float candPtDr = candPt * candDr;
244  size_t icone = std::lower_bound(&cones[0], &cones[ncones], candDr) - &cones[0];
245  float weight2 = candPt * candPt;
246 
247  if (lLead.isNull() || candPt > lLead->pt()) {
248  lSecond = lLead;
249  lLead = icand;
250  } else if ((lSecond.isNull() || candPt > lSecond->pt()) && (candPt < lLead->pt())) {
251  lSecond = icand;
252  }
253 
254  //internalId_.dRMean_ += candPtDr;
255  internalId_.dR2Mean_ += candPtDr * candPtDr;
256 
257  internalId_.ptD_ += candPt * candPt;
258  sumPt += candPt;
259  sumPt2 += candPt * candPt;
260  sum_deta += candDeta * weight2;
261  sum_dphi += candDphi * weight2;
262  sum_deta2 += candDeta * candDeta * weight2;
263  sum_detadphi += candDeta * candDphi * weight2;
264  sum_dphi2 += candDphi * candDphi * weight2;
265  //Teta += candPt * candDR * candDeta;
266  //Tphi += candPt * candDR * candDphi;
267 
268  frac.push_back(candPtFrac);
269  if (icone < ncones) {
270  *coneFracs[icone] += candPt;
271  }
272 
273  if (icand->particleId() == reco::PFCandidate::h0) {
274  if (lLeadNeut.isNull() || candPt > lLeadNeut->pt()) {
275  lLeadNeut = icand;
276  }
277  internalId_.dRMeanNeut_ += candPtDr;
278  fracNeut.push_back(candPtFrac);
279  sumPtNe += candPt;
280  }
281 
282  if (icand->particleId() == reco::PFCandidate::gamma) {
283  if (lLeadEm.isNull() || candPt > lLeadEm->pt()) {
284  lLeadEm = icand;
285  }
286  internalId_.dRMeanEm_ += candPtDr;
287  fracEm.push_back(candPtFrac);
288  sumPtNe += candPt;
289  }
290 
291  if (icand->trackRef().isNonnull() && icand->trackRef().isAvailable()) {
292  if (lLeadCh.isNull() || candPt > lLeadCh->pt()) {
293  lLeadCh = icand;
294  }
295  //internalId_.jetRchg_ += candPtDr;
296  fracCh.push_back(candPtFrac);
297  sumPtCh += candPt;
298  }
299 
300  if (icand->trackRef().isNonnull() && icand->trackRef().isAvailable()) {
301  float tkpt = icand->trackRef()->pt();
302  sumTkPt += tkpt;
303  bool inVtx0 =
304  find(vtx->tracks_begin(), vtx->tracks_end(), reco::TrackBaseRef(icand->trackRef())) != vtx->tracks_end();
305  bool inAnyOther = false;
306 
307  double dZ0 = fabs(icand->trackRef()->dz(vtx->position()));
308  double dZ = dZ0;
309  for (reco::VertexCollection::const_iterator vi = allvtx.begin(); vi != allvtx.end(); ++vi) {
310  const reco::Vertex &iv = *vi;
311  if (iv.isFake() || iv.ndof() < 4) {
312  continue;
313  }
314 
315  bool isVtx0 = (iv.position() - vtx->position()).r() < 0.02;
316 
317  if (!isVtx0 && !inAnyOther) {
318  inAnyOther =
319  find(iv.tracks_begin(), iv.tracks_end(), reco::TrackBaseRef(icand->trackRef())) != iv.tracks_end();
320  }
321 
322  dZ = std::min(dZ, fabs(icand->trackRef()->dz(iv.position())));
323  }
324  if (inVtx0 && !inAnyOther) {
325  internalId_.betaClassic_ += tkpt;
326  } else if (!inVtx0 && inAnyOther) {
327  internalId_.betaStarClassic_ += tkpt;
328  }
329 
330  if (dZ0 < 0.2) {
331  internalId_.beta_ += tkpt;
332  } else if (dZ < 0.2) {
333  internalId_.betaStar_ += tkpt;
334  }
335  }
336 
337  if (lTrail.isNull() || candPt < lTrail->pt()) {
338  lTrail = icand;
339  }
340  }
341 
342  assert(lLead.isNonnull());
343  if (lSecond.isNull()) {
344  lSecond = lTrail;
345  }
346  if (lLeadNeut.isNull()) {
347  lLeadNeut = lTrail;
348  }
349  if (lLeadEm.isNull()) {
350  lLeadEm = lTrail;
351  }
352  if (lLeadCh.isNull()) {
353  lLeadCh = lTrail;
354  }
355  impactTrack = lLeadCh->trackRef();
356 
357  internalId_.nCharged_ = pfjet->chargedMultiplicity();
358  internalId_.nNeutrals_ = pfjet->neutralMultiplicity();
359  internalId_.chgEMfrac_ = pfjet->chargedEmEnergy() / jet->energy();
360  internalId_.neuEMfrac_ = pfjet->neutralEmEnergy() / jet->energy();
361  internalId_.chgHadrfrac_ = pfjet->chargedHadronEnergy() / jet->energy();
362  internalId_.neuHadrfrac_ = pfjet->neutralHadronEnergy() / jet->energy();
363 
364  if (impactTrack.isNonnull() && impactTrack.isAvailable()) {
365  internalId_.d0_ = fabs(impactTrack->dxy(vtx->position()));
366  internalId_.dZ_ = fabs(impactTrack->dz(vtx->position()));
367  } else {
368  internalId_.nParticles_ = constituents.size();
374  std::sort(frac.begin(), frac.end(), std::greater<float>());
375  std::sort(fracCh.begin(), fracCh.end(), std::greater<float>());
376  std::sort(fracEm.begin(), fracEm.end(), std::greater<float>());
377  std::sort(fracNeut.begin(), fracNeut.end(), std::greater<float>());
379 
380  //covMatrix(0,0) /= sumPt2;
381  //covMatrix(0,1) /= sumPt2;
382  //covMatrix(1,1) /= sumPt2;
383  //covMatrix(1,0) = covMatrix(0,1);
384  //internalId_.etaW_ = sqrt(covMatrix(0,0));
385  //internalId_.phiW_ = sqrt(covMatrix(1,1));
386  //internalId_.jetW_ = 0.5*(internalId_.etaW_+internalId_.phiW_);
387  //TVectorD eigVals(2); eigVals = TMatrixDSymEigen(covMatrix).GetEigenValues();
388  //
389  if (internalId_.majW_ < internalId_.minW_) {
390  std::swap(internalId_.majW_, internalId_.minW_);
391  }
392 
393  //internalId_.dRLeadCent_ = reco::deltaR(*jet,*lLead);
394  if (lSecond.isNonnull()) {
395  internalId_.dRLead2nd_ = reco::deltaR(*jet, *lSecond);
396  }
399  //internalId_.jetRchg_ /= jetPt;
400  internalId_.dR2Mean_ /= sumPt2;
401  for (size_t ic = 0; ic < ncones; ++ic) {
402  *coneFracs[ic] /= jetPt;
403  }
404 
405  double ptMean = sumPt / internalId_.nParticles_;
406  double ptRMS = 0;
407  for (unsigned int i0 = 0; i0 < frac.size(); i0++) {
408  ptRMS += (frac[i0] - ptMean) * (frac[i0] - ptMean);
409  }
410  ptRMS /= internalId_.nParticles_;
411  ptRMS = sqrt(ptRMS);
412  internalId_.jetRchg_ = internalId_.leadChPt_ / sumPt;
414 
415  internalId_.ptMean_ = ptMean;
416  internalId_.ptRMS_ = ptRMS / jetPt;
417  internalId_.pt2A_ = sqrt(internalId_.ptD_ / internalId_.nParticles_) / jetPt;
418  internalId_.ptD_ = sqrt(internalId_.ptD_) / sumPt;
420  internalId_.sumChPt_ = sumPtCh;
421  internalId_.sumNePt_ = sumPtNe;
422  if (sumPt > 0) {
423  internalId_.beta_ /= sumPt;
424  internalId_.betaStar_ /= sumPt;
425  } else {
426  assert(internalId_.beta_ == 0. && internalId_.betaStar_ == 0.);
427  }
428  float ave_deta = sum_deta / sumPt2;
429  float ave_dphi = sum_dphi / sumPt2;
430  float ave_deta2 = sum_deta2 / sumPt2;
431  float ave_dphi2 = sum_dphi2 / sumPt2;
432  float a = ave_deta2 - ave_deta * ave_deta;
433  float b = ave_dphi2 - ave_dphi * ave_dphi;
434  float c = -(sum_detadphi / sumPt2 - ave_deta * ave_dphi);
435  float axis1 = 0;
436  float axis2 = 0;
437  if ((((a - b) * (a - b) + 4 * c * c)) > 0) {
438  float delta = sqrt(((a - b) * (a - b) + 4 * c * c));
439  if (a + b + delta > 0) {
440  axis1 = sqrt(0.5 * (a + b + delta));
441  }
442  if (a + b - delta > 0) {
443  axis2 = sqrt(0.5 * (a + b - delta));
444  }
445  } else {
446  axis1 = -1;
447  axis2 = -1;
448  }
449  internalId_.majW_ = axis1; //sqrt(fabs(eigVals(0)));
450  internalId_.minW_ = axis2; //sqrt(fabs(eigVals(1)));
451  //compute Pull
452 
453  float ddetaR_sum(0.0), ddphiR_sum(0.0);
454  for (int i = 0; i < internalId_.nParticles_; ++i) {
456  float weight = part->pt() * part->pt();
457  float deta = part->eta() - jet->eta();
458  float dphi = reco::deltaPhi(*part, *jet);
459  float ddeta, ddphi, ddR;
460  ddeta = deta - ave_deta;
461  ddphi = reco::deltaPhi(dphi, ave_dphi); //2*atan(tan((dphi - ave_dphi)/2.)) ;
462  ddR = sqrt(ddeta * ddeta + ddphi * ddphi);
463  ddetaR_sum += ddR * ddeta * weight;
464  ddphiR_sum += ddR * ddphi * weight;
465  }
466  if (sumPt2 > 0) {
467  float ddetaR_ave = ddetaR_sum / sumPt2;
468  float ddphiR_ave = ddphiR_sum / sumPt2;
469  internalId_.dRMean_ = sqrt(ddetaR_ave * ddetaR_ave + ddphiR_ave * ddphiR_ave);
470  }
471  }
472 
473  if (calculateMva) {
474  runMva();
475  }
476 
478 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
int32_t *__restrict__ iv
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:99
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:152
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:139
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:147
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
Definition: weight.py:1
int neutralMultiplicity() const
neutralMultiplicity
Definition: PFJet.h:154
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
Jets made from PFObjects.
Definition: PFJet.h:20
bool isAvailable() const
Definition: Ptr.h:230
void runMva()
Definition: MVAJetPuId.cc:145
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:41
bool isNull() const
Checks for null.
Definition: Ptr.h:142
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:35
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:19
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
virtual reco::PFCandidatePtr getPFConstituent(unsigned fIndex) const
get specific constituent
Definition: PFJet.cc:27
bool isAvailable() const
Definition: Ref.h:537
void SetPtEtaPhi(const reco::Candidate &p, float &pt, float &eta, float &phi)
Definition: MVAJetPuId.cc:119
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
part
Definition: HCALResponse.h:20
double b
Definition: hdecay.h:118
PileupJetIdentifier internalId_
Definition: MVAJetPuId.h:61
void resetVariables()
Definition: MVAJetPuId.cc:489
double a
Definition: hdecay.h:119
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:95
void Assign(const std::vector< float > &vec, float &a, float &b, float &c, float &d)
Definition: MVAJetPuId.cc:111

◆ computeMva()

PileupJetIdentifier MVAJetPuId::computeMva ( )

Definition at line 191 of file MVAJetPuId.cc.

References internalId_, and runMva().

Referenced by MVAJetPuIdProducer::produce().

191  {
192  runMva();
194 }
void runMva()
Definition: MVAJetPuId.cc:145
PileupJetIdentifier internalId_
Definition: MVAJetPuId.h:61

◆ dumpVariables()

std::string MVAJetPuId::dumpVariables ( ) const

Definition at line 480 of file MVAJetPuId.cc.

References MillePedeFileConverter_cfg::out, and variables_.

480  {
481  std::stringstream out;
482  for (variables_list_t::const_iterator it = variables_.begin(); it != variables_.end(); ++it) {
483  out << std::setw(15) << it->first << std::setw(3) << "=" << std::setw(5) << *it->second.first << " ("
484  << std::setw(5) << it->second.second << ")" << std::endl;
485  }
486  return out.str();
487 }
variables_list_t variables_
Definition: MVAJetPuId.h:62

◆ getJetIdKey()

std::pair< int, int > MVAJetPuId::getJetIdKey ( float  jetPt,
float  jetEta 
)

Definition at line 156 of file MVAJetPuId.cc.

References reco::btau::jetEta, and reco::btau::jetPt.

Referenced by computeIDflag().

156  {
157  int ptId = 0;
158  if (jetPt > 10 && jetPt < 20)
159  ptId = 1;
160  if (jetPt >= 20 && jetPt < 30)
161  ptId = 2;
162  if (jetPt >= 30)
163  ptId = 3;
164 
165  int etaId = 0;
166  if (fabs(jetEta) > 2.5 && fabs(jetEta) < 2.75)
167  etaId = 1;
168  if (fabs(jetEta) >= 2.75 && fabs(jetEta) < 3.0)
169  etaId = 2;
170  if (fabs(jetEta) >= 3.0 && fabs(jetEta) < 5.0)
171  etaId = 3;
172  return std::pair<int, int>(ptId, etaId);
173 }

◆ getVariables()

const variables_list_t& MVAJetPuId::getVariables ( ) const
inline

Definition at line 52 of file MVAJetPuId.h.

References variables_.

52 { return variables_; };
variables_list_t variables_
Definition: MVAJetPuId.h:62

◆ initVariables()

void MVAJetPuId::initVariables ( )
protected

Definition at line 500 of file MVAJetPuId.cc.

References HLT_2022v12_cff::beta, BeamSpotFakeParameters_cfi::betaStar, d0, PileupJetIdentifier::idFlag_, INIT_VARIABLE, internalId_, reco::btau::jetEta, reco::btau::jetPhi, reco::btau::jetPt, LHEJetFilter_cfi::jetR, large_val, beam_dqm_sourceclient-live_cfg::mva, HLT_2022v12_cff::nCharged, jetsAK4_CHS_cff::ptD, rho, and TtFullHadEvtBuilder_cfi::sumPt.

Referenced by setup().

500  {
501  internalId_.idFlag_ = 0;
502  INIT_VARIABLE(mva, "", -100.);
503 
504  INIT_VARIABLE(jetPt, "jetPt", 0.);
505  INIT_VARIABLE(jetEta, "jetEta", large_val);
507  INIT_VARIABLE(jetM, "", 0.);
508  INIT_VARIABLE(nCharged, "nCharged", 0.);
509  INIT_VARIABLE(nNeutrals, "", 0.);
510 
511  INIT_VARIABLE(chgEMfrac, "", 0.);
512  INIT_VARIABLE(neuEMfrac, "", 0.);
513  INIT_VARIABLE(chgHadrfrac, "", 0.);
514  INIT_VARIABLE(neuHadrfrac, "", 0.);
515 
516  INIT_VARIABLE(d0, "", -1000.);
517  INIT_VARIABLE(dZ, "", -1000.);
518  INIT_VARIABLE(nParticles, "nParticles", 0.);
519 
520  INIT_VARIABLE(leadPt, "", 0.);
521  INIT_VARIABLE(leadEta, "", large_val);
522  INIT_VARIABLE(leadPhi, "", large_val);
523  INIT_VARIABLE(secondPt, "", 0.);
524  INIT_VARIABLE(secondEta, "", large_val);
525  INIT_VARIABLE(secondPhi, "", large_val);
526  INIT_VARIABLE(leadNeutPt, "", 0.);
527  INIT_VARIABLE(leadNeutEta, "", large_val);
528 
529  INIT_VARIABLE(jetR, "jetR", 0.);
530  INIT_VARIABLE(pull, "pull", 0.);
531  INIT_VARIABLE(jetRchg, "jetRchg", 0.);
532  INIT_VARIABLE(dR2Mean, "dR2Mean", 0.);
533 
534  INIT_VARIABLE(ptD, "ptD", 0.);
535  INIT_VARIABLE(ptMean, "", 0.);
536  INIT_VARIABLE(ptRMS, "", 0.);
537  INIT_VARIABLE(pt2A, "", 0.);
538  INIT_VARIABLE(ptDCh, "", 0.);
539  INIT_VARIABLE(ptDNe, "", 0.);
540  INIT_VARIABLE(sumPt, "", 0.);
541  INIT_VARIABLE(sumChPt, "", 0.);
542  INIT_VARIABLE(sumNePt, "", 0.);
543  INIT_VARIABLE(secondFrac, "", 0.);
544  INIT_VARIABLE(thirdFrac, "", 0.);
545  INIT_VARIABLE(fourthFrac, "", 0.);
546  INIT_VARIABLE(leadChFrac, "", 0.);
547  INIT_VARIABLE(secondChFrac, "", 0.);
548  INIT_VARIABLE(thirdChFrac, "", 0.);
549  INIT_VARIABLE(fourthChFrac, "", 0.);
550  INIT_VARIABLE(leadNeutFrac, "", 0.);
551  INIT_VARIABLE(secondNeutFrac, "", 0.);
552  INIT_VARIABLE(thirdNeutFrac, "", 0.);
553  INIT_VARIABLE(fourthNeutFrac, "", 0.);
554  INIT_VARIABLE(leadEmFrac, "", 0.);
555  INIT_VARIABLE(secondEmFrac, "", 0.);
556  INIT_VARIABLE(thirdEmFrac, "", 0.);
557  INIT_VARIABLE(fourthEmFrac, "", 0.);
558  INIT_VARIABLE(jetW, "", 1.);
559  INIT_VARIABLE(etaW, "", 1.);
560  INIT_VARIABLE(phiW, "", 1.);
561  INIT_VARIABLE(majW, "majW", 1.);
562  INIT_VARIABLE(minW, "minW", 1.);
563  INIT_VARIABLE(frac01, "frac01", 0.);
564  INIT_VARIABLE(frac02, "frac02", 0.);
565  INIT_VARIABLE(frac03, "frac03", 0.);
566  INIT_VARIABLE(frac04, "frac04", 0.);
567 
568  INIT_VARIABLE(beta, "beta", 0.);
569  INIT_VARIABLE(betaStar, "betaStar", 0.);
570  INIT_VARIABLE(betaClassic, "betaClassic", 0.);
571  INIT_VARIABLE(betaStarClassic, "betaStarClassic", 0.);
572  INIT_VARIABLE(rho, "rho", 0.);
573 }
#define INIT_VARIABLE(NAME, TMVANAME, VAL)
Definition: MVAJetPuId.cc:496
static constexpr float d0
const float large_val
Definition: MVAJetPuId.cc:14
PileupJetIdentifier internalId_
Definition: MVAJetPuId.h:61

◆ method()

const std::string MVAJetPuId::method ( ) const
inline

Definition at line 41 of file MVAJetPuId.h.

References tmvaMethod_.

41 { return tmvaMethod_; }
std::string tmvaMethod_
Definition: MVAJetPuId.h:65

◆ resetVariables()

void MVAJetPuId::resetVariables ( )
protected

Definition at line 489 of file MVAJetPuId.cc.

References PileupJetIdentifier::idFlag_, internalId_, and variables_.

Referenced by computeIdVariables().

489  {
490  internalId_.idFlag_ = 0;
491  for (variables_list_t::iterator it = variables_.begin(); it != variables_.end(); ++it) {
492  *it->second.first = it->second.second;
493  }
494 }
variables_list_t variables_
Definition: MVAJetPuId.h:62
PileupJetIdentifier internalId_
Definition: MVAJetPuId.h:61

◆ runMva()

void MVAJetPuId::runMva ( )
protected

Definition at line 145 of file MVAJetPuId.cc.

References bookReader(), computeIDflag(), PileupJetIdentifier::idFlag_, internalId_, PileupJetIdentifier::mva_, reader_, and tmvaMethod_.

Referenced by computeIdVariables(), and computeMva().

145  {
146  if (!reader_) {
147  bookReader();
148  }
149  if (fabs(internalId_.jetEta_) < 5.0)
150  internalId_.mva_ = reader_->EvaluateMVA(tmvaMethod_.c_str());
151  if (fabs(internalId_.jetEta_) >= 5.0)
152  internalId_.mva_ = -2.;
154 }
void bookReader()
Definition: MVAJetPuId.cc:125
TMVA::Reader * reader_
Definition: MVAJetPuId.h:64
PileupJetIdentifier internalId_
Definition: MVAJetPuId.h:61
int computeIDflag(float mva, float jetPt, float jetEta)
Definition: MVAJetPuId.cc:175
std::string tmvaMethod_
Definition: MVAJetPuId.h:65

◆ set()

void MVAJetPuId::set ( const PileupJetIdentifier id)

Definition at line 143 of file MVAJetPuId.cc.

References triggerObjects_cff::id, and internalId_.

Referenced by MVAJetPuIdProducer::produce().

143 { internalId_ = id; }
PileupJetIdentifier internalId_
Definition: MVAJetPuId.h:61

◆ setup()

void MVAJetPuId::setup ( )
protected

Definition at line 66 of file MVAJetPuId.cc.

References initVariables(), tmvaNames_, and tmvaVariables_.

Referenced by MVAJetPuId().

66  {
67  initVariables();
68 
69  tmvaVariables_.clear();
70  tmvaVariables_.push_back("rho");
71  tmvaVariables_.push_back("nParticles");
72  tmvaVariables_.push_back("nCharged");
73  tmvaVariables_.push_back("majW");
74  tmvaVariables_.push_back("minW");
75  tmvaVariables_.push_back("frac01");
76  tmvaVariables_.push_back("frac02");
77  tmvaVariables_.push_back("frac03");
78  tmvaVariables_.push_back("frac04");
79  tmvaVariables_.push_back("ptD");
80  tmvaVariables_.push_back("beta");
81  tmvaVariables_.push_back("betaStar");
82  tmvaVariables_.push_back("dR2Mean");
83  tmvaVariables_.push_back("pull");
84  tmvaVariables_.push_back("jetR");
85  tmvaVariables_.push_back("jetRchg");
86 
87  tmvaNames_["rho"] = "rho";
88  tmvaNames_["nParticles"] = "nParticles";
89  tmvaNames_["nCharged"] = "nCharged";
90  tmvaNames_["majW"] = "majW";
91  tmvaNames_["minW"] = "minW";
92  tmvaNames_["frac01"] = "frac01";
93  tmvaNames_["frac02"] = "frac02";
94  tmvaNames_["frac03"] = "frac03";
95  tmvaNames_["frac04"] = "frac04";
96  tmvaNames_["ptD"] = "ptD";
97  tmvaNames_["beta"] = "beta";
98  tmvaNames_["betaStar"] = "betaStar";
99  tmvaNames_["dR2Mean"] = "dR2Mean";
100  tmvaNames_["pull"] = "pull";
101  tmvaNames_["jetR"] = "jetR";
102  tmvaNames_["jetRchg"] = "jetRchg";
103 }
std::vector< std::string > tmvaVariables_
Definition: MVAJetPuId.h:66
void initVariables()
Definition: MVAJetPuId.cc:500
std::map< std::string, std::string > tmvaNames_
Definition: MVAJetPuId.h:68

Member Data Documentation

◆ betaStarCut_

Float_t MVAJetPuId::betaStarCut_[NWPs][NEtas][NPts]
protected

Definition at line 75 of file MVAJetPuId.h.

◆ cutBased_

bool MVAJetPuId::cutBased_
protected

Definition at line 72 of file MVAJetPuId.h.

◆ impactParTkThreshod_

Float_t MVAJetPuId::impactParTkThreshod_
protected

Definition at line 71 of file MVAJetPuId.h.

Referenced by MVAJetPuId().

◆ internalId_

PileupJetIdentifier MVAJetPuId::internalId_
protected

Definition at line 61 of file MVAJetPuId.h.

Referenced by computeIdVariables(), computeMva(), initVariables(), resetVariables(), runMva(), and set().

◆ mvacut_

Float_t MVAJetPuId::mvacut_[NWPs][NEtas][NPts]
protected

Definition at line 73 of file MVAJetPuId.h.

Referenced by computeIDflag(), and MVAJetPuId().

◆ NEtas

constexpr int MVAJetPuId::NEtas = 4
staticprivate

Definition at line 19 of file MVAJetPuId.h.

◆ NPts

constexpr int MVAJetPuId::NPts = 4
staticprivate

Definition at line 18 of file MVAJetPuId.h.

Referenced by MVAJetPuId().

◆ NWPs

constexpr int MVAJetPuId::NWPs = 3
staticprivate

Definition at line 17 of file MVAJetPuId.h.

Referenced by MVAJetPuId().

◆ reader_

TMVA::Reader* MVAJetPuId::reader_
protected

Definition at line 64 of file MVAJetPuId.h.

Referenced by bookReader(), MVAJetPuId(), runMva(), and ~MVAJetPuId().

◆ rmsCut_

Float_t MVAJetPuId::rmsCut_[NWPs][NEtas][NPts]
protected

Definition at line 74 of file MVAJetPuId.h.

◆ tmvaMethod_

std::string MVAJetPuId::tmvaMethod_
protected

Definition at line 65 of file MVAJetPuId.h.

Referenced by bookReader(), method(), MVAJetPuId(), and runMva().

◆ tmvaNames_

std::map<std::string, std::string> MVAJetPuId::tmvaNames_
protected

Definition at line 68 of file MVAJetPuId.h.

Referenced by bookReader(), and setup().

◆ tmvaSpectators_

std::vector<std::string> MVAJetPuId::tmvaSpectators_
protected

Definition at line 67 of file MVAJetPuId.h.

Referenced by bookReader(), and MVAJetPuId().

◆ tmvaVariables_

std::vector<std::string> MVAJetPuId::tmvaVariables_
protected

Definition at line 66 of file MVAJetPuId.h.

Referenced by bookReader(), MVAJetPuId(), and setup().

◆ tmvaWeights_

std::string MVAJetPuId::tmvaWeights_
protected

Definition at line 65 of file MVAJetPuId.h.

Referenced by bookReader(), and MVAJetPuId().

◆ variables_

variables_list_t MVAJetPuId::variables_
protected

Definition at line 62 of file MVAJetPuId.h.

Referenced by bookReader(), dumpVariables(), getVariables(), and resetVariables().

◆ version_

Int_t MVAJetPuId::version_
protected

Definition at line 70 of file MVAJetPuId.h.

Referenced by MVAJetPuId().