Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-08-05 22:48:31

0001 #include "RecoJets/JetProducers/interface/PileupJetIdAlgo.h"
0002 
0003 #include "DataFormats/JetReco/interface/PFJet.h"
0004 #include "DataFormats/JetReco/interface/Jet.h"
0005 #include "DataFormats/VertexReco/interface/Vertex.h"
0006 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0007 #include "DataFormats/Math/interface/deltaR.h"
0008 #include "FWCore/ParameterSet/interface/FileInPath.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include "CommonTools/MVAUtils/interface/GBRForestTools.h"
0011 
0012 #include "TMatrixDSym.h"
0013 #include "TMatrixDSymEigen.h"
0014 
0015 #include <utility>
0016 
0017 // ------------------------------------------------------------------------------------------
0018 const float large_val = std::numeric_limits<float>::max();
0019 
0020 // ------------------------------------------------------------------------------------------
0021 
0022 PileupJetIdAlgo::AlgoGBRForestsAndConstants::AlgoGBRForestsAndConstants(edm::ParameterSet const& ps, bool runMvas)
0023     : cutBased_(ps.getParameter<bool>("cutBased")),
0024       etaBinnedWeights_(false),
0025       runMvas_(runMvas),
0026       nEtaBins_(0),
0027       label_(ps.getParameter<std::string>("label")),
0028       mvacut_{},
0029       rmsCut_{},
0030       betaStarCut_{} {
0031   std::vector<edm::FileInPath> tmvaEtaWeights;
0032   std::vector<std::string> tmvaSpectators;
0033   int version;
0034 
0035   if (!cutBased_) {
0036     etaBinnedWeights_ = ps.getParameter<bool>("etaBinnedWeights");
0037     if (etaBinnedWeights_) {
0038       const std::vector<edm::ParameterSet>& trainings = ps.getParameter<std::vector<edm::ParameterSet> >("trainings");
0039       nEtaBins_ = ps.getParameter<int>("nEtaBins");
0040       for (int v = 0; v < nEtaBins_; v++) {
0041         tmvaEtaWeights.push_back(trainings.at(v).getParameter<edm::FileInPath>("tmvaWeights"));
0042         jEtaMin_.push_back(trainings.at(v).getParameter<double>("jEtaMin"));
0043         jEtaMax_.push_back(trainings.at(v).getParameter<double>("jEtaMax"));
0044       }
0045       for (int v = 0; v < nEtaBins_; v++) {
0046         tmvaEtaVariables_.push_back(trainings.at(v).getParameter<std::vector<std::string> >("tmvaVariables"));
0047       }
0048     } else {
0049       tmvaVariables_ = ps.getParameter<std::vector<std::string> >("tmvaVariables");
0050     }
0051     tmvaMethod_ = ps.getParameter<std::string>("tmvaMethod");
0052     tmvaSpectators = ps.getParameter<std::vector<std::string> >("tmvaSpectators");
0053     version = ps.getParameter<int>("version");
0054   } else {
0055     version = USER;
0056   }
0057 
0058   edm::ParameterSet jetConfig = ps.getParameter<edm::ParameterSet>("JetIdParams");
0059   for (int i0 = 0; i0 < 3; i0++) {
0060     std::string lCutType = "Tight";
0061     if (i0 == PileupJetIdentifier::kMedium)
0062       lCutType = "Medium";
0063     if (i0 == PileupJetIdentifier::kLoose)
0064       lCutType = "Loose";
0065     int nCut = 1;
0066     if (cutBased_)
0067       nCut++;
0068     for (int i1 = 0; i1 < nCut; i1++) {
0069       std::string lFullCutType = lCutType;
0070       if (cutBased_ && i1 == 0)
0071         lFullCutType = "BetaStar" + lCutType;
0072       if (cutBased_ && i1 == 1)
0073         lFullCutType = "RMS" + lCutType;
0074       std::vector<double> pt010 = jetConfig.getParameter<std::vector<double> >(("Pt010_" + lFullCutType).c_str());
0075       std::vector<double> pt1020 = jetConfig.getParameter<std::vector<double> >(("Pt1020_" + lFullCutType).c_str());
0076       std::vector<double> pt2030 = jetConfig.getParameter<std::vector<double> >(("Pt2030_" + lFullCutType).c_str());
0077       std::vector<double> pt3040 = jetConfig.getParameter<std::vector<double> >(("Pt3040_" + lFullCutType).c_str());
0078       std::vector<double> pt4050 = jetConfig.getParameter<std::vector<double> >(("Pt4050_" + lFullCutType).c_str());
0079       if (!cutBased_) {
0080         for (int i2 = 0; i2 < 4; i2++)
0081           mvacut_[i0][0][i2] = pt010[i2];
0082         for (int i2 = 0; i2 < 4; i2++)
0083           mvacut_[i0][1][i2] = pt1020[i2];
0084         for (int i2 = 0; i2 < 4; i2++)
0085           mvacut_[i0][2][i2] = pt2030[i2];
0086         for (int i2 = 0; i2 < 4; i2++)
0087           mvacut_[i0][3][i2] = pt3040[i2];
0088         for (int i2 = 0; i2 < 4; i2++)
0089           mvacut_[i0][4][i2] = pt4050[i2];
0090       }
0091       if (cutBased_ && i1 == 0) {
0092         for (int i2 = 0; i2 < 4; i2++)
0093           betaStarCut_[i0][0][i2] = pt010[i2];
0094         for (int i2 = 0; i2 < 4; i2++)
0095           betaStarCut_[i0][1][i2] = pt1020[i2];
0096         for (int i2 = 0; i2 < 4; i2++)
0097           betaStarCut_[i0][2][i2] = pt2030[i2];
0098         for (int i2 = 0; i2 < 4; i2++)
0099           betaStarCut_[i0][3][i2] = pt3040[i2];
0100         for (int i2 = 0; i2 < 4; i2++)
0101           betaStarCut_[i0][4][i2] = pt4050[i2];
0102       }
0103       if (cutBased_ && i1 == 1) {
0104         for (int i2 = 0; i2 < 4; i2++)
0105           rmsCut_[i0][0][i2] = pt010[i2];
0106         for (int i2 = 0; i2 < 4; i2++)
0107           rmsCut_[i0][1][i2] = pt1020[i2];
0108         for (int i2 = 0; i2 < 4; i2++)
0109           rmsCut_[i0][2][i2] = pt2030[i2];
0110         for (int i2 = 0; i2 < 4; i2++)
0111           rmsCut_[i0][3][i2] = pt3040[i2];
0112         for (int i2 = 0; i2 < 4; i2++)
0113           rmsCut_[i0][4][i2] = pt4050[i2];
0114       }
0115     }
0116   }
0117 
0118   if (!cutBased_) {
0119     assert(tmvaMethod_.empty() || ((!tmvaVariables_.empty() || (!tmvaEtaVariables_.empty())) && version == USER));
0120   }
0121 
0122   if ((!cutBased_) && (runMvas_)) {
0123     if (etaBinnedWeights_) {
0124       for (int v = 0; v < nEtaBins_; v++) {
0125         etaReader_.push_back(createGBRForest(tmvaEtaWeights.at(v)));
0126       }
0127     } else {
0128       reader_ = createGBRForest(ps.getParameter<std::string>("tmvaWeights"));
0129     }
0130   }
0131 }
0132 
0133 PileupJetIdAlgo::PileupJetIdAlgo(AlgoGBRForestsAndConstants const* cache) : cache_(cache) { initVariables(); }
0134 
0135 // ------------------------------------------------------------------------------------------
0136 PileupJetIdAlgo::~PileupJetIdAlgo() {}
0137 
0138 // ------------------------------------------------------------------------------------------
0139 void assign(const std::vector<float>& vec, float& a, float& b, float& c, float& d) {
0140   size_t sz = vec.size();
0141   a = (sz > 0 ? vec[0] : 0.);
0142   b = (sz > 1 ? vec[1] : 0.);
0143   c = (sz > 2 ? vec[2] : 0.);
0144   d = (sz > 3 ? vec[3] : 0.);
0145 }
0146 // ------------------------------------------------------------------------------------------
0147 void setPtEtaPhi(const reco::Candidate& p, float& pt, float& eta, float& phi) {
0148   pt = p.pt();
0149   eta = p.eta();
0150   phi = p.phi();
0151 }
0152 
0153 // ------------------------------------------------------------------------------------------
0154 void PileupJetIdAlgo::set(const PileupJetIdentifier& id) { internalId_ = id; }
0155 
0156 // ------------------------------------------------------------------------------------------
0157 
0158 float PileupJetIdAlgo::getMVAval(const std::vector<std::string>& varList,
0159                                  const std::unique_ptr<const GBRForest>& reader) {
0160   std::vector<float> vars;
0161   for (std::vector<std::string>::const_iterator it = varList.begin(); it != varList.end(); ++it) {
0162     std::pair<float*, float> var = variables_.at(*it);
0163     vars.push_back(*var.first);
0164   }
0165   return reader->GetClassifier(vars.data());
0166 }
0167 
0168 void PileupJetIdAlgo::runMva() {
0169   if (cache_->cutBased()) {
0170     internalId_.idFlag_ = computeCutIDflag(
0171         internalId_.betaStarClassic_, internalId_.dR2Mean_, internalId_.nvtx_, internalId_.jetPt_, internalId_.jetEta_);
0172   } else {
0173     if (std::abs(internalId_.jetEta_) >= 5.0) {
0174       internalId_.mva_ = -2.;
0175     } else {
0176       if (cache_->etaBinnedWeights()) {
0177         if (std::abs(internalId_.jetEta_) > cache_->jEtaMax().at(cache_->nEtaBins() - 1)) {
0178           internalId_.mva_ = -2.;
0179         } else {
0180           for (int v = 0; v < cache_->nEtaBins(); v++) {
0181             if (std::abs(internalId_.jetEta_) >= cache_->jEtaMin().at(v) &&
0182                 std::abs(internalId_.jetEta_) < cache_->jEtaMax().at(v)) {
0183               internalId_.mva_ = getMVAval(cache_->tmvaEtaVariables().at(v), cache_->etaReader().at(v));
0184               break;
0185             }
0186           }
0187         }
0188       } else {
0189         internalId_.mva_ = getMVAval(cache_->tmvaVariables(), cache_->reader());
0190       }
0191     }
0192     internalId_.idFlag_ = computeIDflag(internalId_.mva_, internalId_.jetPt_, internalId_.jetEta_);
0193   }
0194 }
0195 
0196 // ------------------------------------------------------------------------------------------
0197 std::pair<int, int> PileupJetIdAlgo::getJetIdKey(float jetPt, float jetEta) {
0198   int ptId = 0;
0199   if (jetPt >= 10 && jetPt < 20)
0200     ptId = 1;
0201   if (jetPt >= 20 && jetPt < 30)
0202     ptId = 2;
0203   if (jetPt >= 30 && jetPt < 40)
0204     ptId = 3;
0205   if (jetPt >= 40)
0206     ptId = 4;
0207 
0208   int etaId = 0;
0209   if (std::abs(jetEta) >= 2.5 && std::abs(jetEta) < 2.75)
0210     etaId = 1;
0211   if (std::abs(jetEta) >= 2.75 && std::abs(jetEta) < 3.0)
0212     etaId = 2;
0213   if (std::abs(jetEta) >= 3.0 && std::abs(jetEta) < 5.0)
0214     etaId = 3;
0215 
0216   return std::pair<int, int>(ptId, etaId);
0217 }
0218 // ------------------------------------------------------------------------------------------
0219 int PileupJetIdAlgo::computeCutIDflag(float betaStarClassic, float dR2Mean, float nvtx, float jetPt, float jetEta) {
0220   std::pair<int, int> jetIdKey = getJetIdKey(jetPt, jetEta);
0221   float betaStarModified = betaStarClassic / log(nvtx - 0.64);
0222   int idFlag(0);
0223   if (betaStarModified < cache_->betaStarCut()[PileupJetIdentifier::kTight][jetIdKey.first][jetIdKey.second] &&
0224       dR2Mean < cache_->rmsCut()[PileupJetIdentifier::kTight][jetIdKey.first][jetIdKey.second])
0225     idFlag += 1 << PileupJetIdentifier::kTight;
0226 
0227   if (betaStarModified < cache_->betaStarCut()[PileupJetIdentifier::kMedium][jetIdKey.first][jetIdKey.second] &&
0228       dR2Mean < cache_->rmsCut()[PileupJetIdentifier::kMedium][jetIdKey.first][jetIdKey.second])
0229     idFlag += 1 << PileupJetIdentifier::kMedium;
0230 
0231   if (betaStarModified < cache_->betaStarCut()[PileupJetIdentifier::kLoose][jetIdKey.first][jetIdKey.second] &&
0232       dR2Mean < cache_->rmsCut()[PileupJetIdentifier::kLoose][jetIdKey.first][jetIdKey.second])
0233     idFlag += 1 << PileupJetIdentifier::kLoose;
0234   return idFlag;
0235 }
0236 // ------------------------------------------------------------------------------------------
0237 int PileupJetIdAlgo::computeIDflag(float mva, float jetPt, float jetEta) {
0238   std::pair<int, int> jetIdKey = getJetIdKey(jetPt, jetEta);
0239   return computeIDflag(mva, jetIdKey.first, jetIdKey.second);
0240 }
0241 
0242 // ------------------------------------------------------------------------------------------
0243 int PileupJetIdAlgo::computeIDflag(float mva, int ptId, int etaId) {
0244   int idFlag(0);
0245   if (mva > cache_->mvacut()[PileupJetIdentifier::kTight][ptId][etaId])
0246     idFlag += 1 << PileupJetIdentifier::kTight;
0247   if (mva > cache_->mvacut()[PileupJetIdentifier::kMedium][ptId][etaId])
0248     idFlag += 1 << PileupJetIdentifier::kMedium;
0249   if (mva > cache_->mvacut()[PileupJetIdentifier::kLoose][ptId][etaId])
0250     idFlag += 1 << PileupJetIdentifier::kLoose;
0251   return idFlag;
0252 }
0253 
0254 // ------------------------------------------------------------------------------------------
0255 PileupJetIdentifier PileupJetIdAlgo::computeMva() {
0256   runMva();
0257   return PileupJetIdentifier(internalId_);
0258 }
0259 
0260 // ------------------------------------------------------------------------------------------
0261 PileupJetIdentifier PileupJetIdAlgo::computeIdVariables(const reco::Jet* jet,
0262                                                         float jec,
0263                                                         const reco::Vertex* vtx,
0264                                                         const reco::VertexCollection& allvtx,
0265                                                         double rho,
0266                                                         bool usePuppi) {
0267   // initialize all variables to 0
0268   resetVariables();
0269 
0270   // loop over constituents, accumulate sums and find leading candidates
0271   const pat::Jet* patjet = dynamic_cast<const pat::Jet*>(jet);
0272   const reco::PFJet* pfjet = dynamic_cast<const reco::PFJet*>(jet);
0273   assert(patjet != nullptr || pfjet != nullptr);
0274   if (patjet != nullptr && jec == 0.) {  // if this is a pat jet and no jec has been passed take the jec from the object
0275     jec = patjet->pt() / patjet->correctedJet(0).pt();
0276   }
0277   if (jec <= 0.) {
0278     jec = 1.;
0279   }
0280 
0281   const reco::Candidate *lLead = nullptr, *lSecond = nullptr, *lLeadNeut = nullptr, *lLeadEm = nullptr,
0282                         *lLeadCh = nullptr, *lTrail = nullptr;
0283   std::vector<float> frac, fracCh, fracEm, fracNeut;
0284   float cones[] = {0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7};
0285   size_t ncones = sizeof(cones) / sizeof(float);
0286   float* coneFracs[] = {&internalId_.frac01_,
0287                         &internalId_.frac02_,
0288                         &internalId_.frac03_,
0289                         &internalId_.frac04_,
0290                         &internalId_.frac05_,
0291                         &internalId_.frac06_,
0292                         &internalId_.frac07_};
0293   float* coneEmFracs[] = {&internalId_.emFrac01_,
0294                           &internalId_.emFrac02_,
0295                           &internalId_.emFrac03_,
0296                           &internalId_.emFrac04_,
0297                           &internalId_.emFrac05_,
0298                           &internalId_.emFrac06_,
0299                           &internalId_.emFrac07_};
0300   float* coneNeutFracs[] = {&internalId_.neutFrac01_,
0301                             &internalId_.neutFrac02_,
0302                             &internalId_.neutFrac03_,
0303                             &internalId_.neutFrac04_,
0304                             &internalId_.neutFrac05_,
0305                             &internalId_.neutFrac06_,
0306                             &internalId_.neutFrac07_};
0307   float* coneChFracs[] = {&internalId_.chFrac01_,
0308                           &internalId_.chFrac02_,
0309                           &internalId_.chFrac03_,
0310                           &internalId_.chFrac04_,
0311                           &internalId_.chFrac05_,
0312                           &internalId_.chFrac06_,
0313                           &internalId_.chFrac07_};
0314   TMatrixDSym covMatrix(2);
0315   covMatrix = 0.;
0316   float jetPt = jet->pt() / jec;  // use uncorrected pt for shape variables
0317   float sumPt = 0., sumPt2 = 0., sumTkPt = 0., sumPtCh = 0, sumPtNe = 0;
0318   float multNeut = 0.0;
0319   setPtEtaPhi(
0320       *jet, internalId_.jetPt_, internalId_.jetEta_, internalId_.jetPhi_);  // use corrected pt for jet kinematics
0321   internalId_.jetM_ = jet->mass();
0322   internalId_.nvtx_ = allvtx.size();
0323   internalId_.rho_ = rho;
0324 
0325   float dRmin(1000);
0326 
0327   for (unsigned i = 0; i < jet->numberOfSourceCandidatePtrs(); ++i) {
0328     reco::CandidatePtr pfJetConstituent = jet->sourceCandidatePtr(i);
0329 
0330     const reco::Candidate* icand = pfJetConstituent.get();
0331     const pat::PackedCandidate* lPack = dynamic_cast<const pat::PackedCandidate*>(icand);
0332     const reco::PFCandidate* lPF = dynamic_cast<const reco::PFCandidate*>(icand);
0333     bool isPacked = true;
0334     if (lPack == nullptr) {
0335       isPacked = false;
0336     }
0337     float candPuppiWeight = 1.0;
0338     if (usePuppi && isPacked)
0339       candPuppiWeight = lPack->puppiWeight();
0340     float candPt = (icand->pt()) * candPuppiWeight;
0341     float candPtFrac = candPt / jetPt;
0342     float candDr = reco::deltaR(*icand, *jet);
0343     float candDeta = icand->eta() - jet->eta();
0344     float candDphi = reco::deltaPhi(*icand, *jet);
0345     float candPtDr = candPt * candDr;
0346     size_t icone = std::lower_bound(&cones[0], &cones[ncones], candDr) - &cones[0];
0347 
0348     if (candDr < dRmin)
0349       dRmin = candDr;
0350 
0351     // // all particles
0352     if (lLead == nullptr || candPt > lLead->pt()) {
0353       lSecond = lLead;
0354       lLead = icand;
0355     } else if ((lSecond == nullptr || candPt > lSecond->pt()) && (candPt < lLead->pt())) {
0356       lSecond = icand;
0357     }
0358 
0359     // // average shapes
0360     internalId_.dRMean_ += candPtDr;
0361     internalId_.dR2Mean_ += candPtDr * candPtDr;
0362     covMatrix(0, 0) += candPt * candPt * candDeta * candDeta;
0363     covMatrix(0, 1) += candPt * candPt * candDeta * candDphi;
0364     covMatrix(1, 1) += candPt * candPt * candDphi * candDphi;
0365     internalId_.ptD_ += candPt * candPt;
0366     sumPt += candPt;
0367     sumPt2 += candPt * candPt;
0368 
0369     // single most energetic candiates and jet shape profiles
0370     frac.push_back(candPtFrac);
0371 
0372     if (icone < ncones) {
0373       *coneFracs[icone] += candPt;
0374     }
0375 
0376     // neutrals
0377     if (abs(icand->pdgId()) == 130) {
0378       if (lLeadNeut == nullptr || candPt > lLeadNeut->pt()) {
0379         lLeadNeut = icand;
0380       }
0381       internalId_.dRMeanNeut_ += candPtDr;
0382       fracNeut.push_back(candPtFrac);
0383       if (icone < ncones) {
0384         *coneNeutFracs[icone] += candPt;
0385       }
0386       internalId_.ptDNe_ += candPt * candPt;
0387       sumPtNe += candPt;
0388       multNeut += candPuppiWeight;
0389     }
0390     // EM candidated
0391     if (icand->pdgId() == 22) {
0392       if (lLeadEm == nullptr || candPt > lLeadEm->pt()) {
0393         lLeadEm = icand;
0394       }
0395       internalId_.dRMeanEm_ += candPtDr;
0396       fracEm.push_back(candPtFrac);
0397       if (icone < ncones) {
0398         *coneEmFracs[icone] += candPt;
0399       }
0400       internalId_.ptDNe_ += candPt * candPt;
0401       sumPtNe += candPt;
0402       multNeut += candPuppiWeight;
0403     }
0404     if ((abs(icand->pdgId()) == 1) || (abs(icand->pdgId()) == 2))
0405       multNeut += candPuppiWeight;
0406 
0407     // Charged  particles
0408     if (icand->charge() != 0) {
0409       if (lLeadCh == nullptr || candPt > lLeadCh->pt()) {
0410         lLeadCh = icand;
0411 
0412         const reco::Track* pfTrk = icand->bestTrack();
0413         if (lPF && std::abs(icand->pdgId()) == 13 && pfTrk == nullptr) {
0414           reco::MuonRef lmuRef = lPF->muonRef();
0415           if (lmuRef.isNonnull()) {
0416             const reco::Muon& lmu = *lmuRef.get();
0417             pfTrk = lmu.bestTrack();
0418             edm::LogWarning("BadMuon")
0419                 << "Found a PFCandidate muon without a trackRef: falling back to Muon::bestTrack ";
0420           }
0421         }
0422         if (pfTrk == nullptr) {  //protection against empty pointers for the miniAOD case
0423           //To handle the electron case
0424           if (isPacked) {
0425             internalId_.d0_ = std::abs(lPack->dxy(vtx->position()));
0426             internalId_.dZ_ = std::abs(lPack->dz(vtx->position()));
0427           } else if (lPF != nullptr) {
0428             pfTrk = (lPF->trackRef().get() == nullptr) ? lPF->gsfTrackRef().get() : lPF->trackRef().get();
0429             internalId_.d0_ = std::abs(pfTrk->dxy(vtx->position()));
0430             internalId_.dZ_ = std::abs(pfTrk->dz(vtx->position()));
0431           }
0432         } else {
0433           internalId_.d0_ = std::abs(pfTrk->dxy(vtx->position()));
0434           internalId_.dZ_ = std::abs(pfTrk->dz(vtx->position()));
0435         }
0436       }
0437       internalId_.dRMeanCh_ += candPtDr;
0438       internalId_.ptDCh_ += candPt * candPt;
0439       fracCh.push_back(candPtFrac);
0440       if (icone < ncones) {
0441         *coneChFracs[icone] += candPt;
0442       }
0443       sumPtCh += candPt;
0444     }
0445     // // beta and betastar
0446     if (icand->charge() != 0) {
0447       if (!isPacked) {
0448         if (lPF->trackRef().isNonnull()) {
0449           float tkpt = candPt;
0450           sumTkPt += tkpt;
0451           // 'classic' beta definition based on track-vertex association
0452           bool inVtx0 = vtx->trackWeight(lPF->trackRef()) > 0;
0453 
0454           bool inAnyOther = false;
0455           // alternative beta definition based on track-vertex distance of closest approach
0456           double dZ0 = std::abs(lPF->trackRef()->dz(vtx->position()));
0457           double dZ = dZ0;
0458           for (reco::VertexCollection::const_iterator vi = allvtx.begin(); vi != allvtx.end(); ++vi) {
0459             const reco::Vertex& iv = *vi;
0460             if (iv.isFake() || iv.ndof() < 4) {
0461               continue;
0462             }
0463             // the primary vertex may have been copied by the user: check identity by position
0464             bool isVtx0 = (iv.position() - vtx->position()).r() < 0.02;
0465             // 'classic' beta definition: check if the track is associated with
0466             // any vertex other than the primary one
0467             if (!isVtx0 && !inAnyOther) {
0468               inAnyOther = vtx->trackWeight(lPF->trackRef()) <= 0;
0469             }
0470             // alternative beta: find closest vertex to the track
0471             dZ = std::min(dZ, std::abs(lPF->trackRef()->dz(iv.position())));
0472           }
0473           // classic beta/betaStar
0474           if (inVtx0 && !inAnyOther) {
0475             internalId_.betaClassic_ += tkpt;
0476           } else if (!inVtx0 && inAnyOther) {
0477             internalId_.betaStarClassic_ += tkpt;
0478           }
0479           // alternative beta/betaStar
0480           if (dZ0 < 0.2) {
0481             internalId_.beta_ += tkpt;
0482           } else if (dZ < 0.2) {
0483             internalId_.betaStar_ += tkpt;
0484           }
0485         }
0486       } else {
0487         float tkpt = candPt;
0488         sumTkPt += tkpt;
0489         bool inVtx0 = false;
0490         bool inVtxOther = false;
0491         double dZ0 = 9999.;
0492         double dZ_tmp = 9999.;
0493         for (unsigned vtx_i = 0; vtx_i < allvtx.size(); vtx_i++) {
0494           auto iv = allvtx[vtx_i];
0495 
0496           if (iv.isFake())
0497             continue;
0498 
0499           // Match to vertex in case of copy as above
0500           bool isVtx0 = (iv.position() - vtx->position()).r() < 0.02;
0501 
0502           if (isVtx0) {
0503             if (lPack->fromPV(vtx_i) == pat::PackedCandidate::PVUsedInFit)
0504               inVtx0 = true;
0505             if (lPack->fromPV(vtx_i) == 0)
0506               inVtxOther = true;
0507             dZ0 = lPack->dz(iv.position());
0508           }
0509 
0510           if (fabs(lPack->dz(iv.position())) < fabs(dZ_tmp)) {
0511             dZ_tmp = lPack->dz(iv.position());
0512           }
0513         }
0514         if (inVtx0) {
0515           internalId_.betaClassic_ += tkpt;
0516         } else if (inVtxOther) {
0517           internalId_.betaStarClassic_ += tkpt;
0518         }
0519         if (fabs(dZ0) < 0.2) {
0520           internalId_.beta_ += tkpt;
0521         } else if (fabs(dZ_tmp) < 0.2) {
0522           internalId_.betaStar_ += tkpt;
0523         }
0524       }
0525     }
0526 
0527     // trailing candidate
0528     if (lTrail == nullptr || candPt < lTrail->pt()) {
0529       lTrail = icand;
0530     }
0531   }
0532 
0533   // // Finalize all variables
0534   assert(!(lLead == nullptr));
0535 
0536   if (lSecond == nullptr) {
0537     lSecond = lTrail;
0538   }
0539   if (lLeadNeut == nullptr) {
0540     lLeadNeut = lTrail;
0541   }
0542   if (lLeadEm == nullptr) {
0543     lLeadEm = lTrail;
0544   }
0545   if (lLeadCh == nullptr) {
0546     lLeadCh = lTrail;
0547   }
0548 
0549   if (patjet != nullptr) {  // to enable running on MiniAOD slimmedJets
0550     internalId_.nCharged_ = patjet->chargedMultiplicity();
0551     internalId_.nNeutrals_ = patjet->neutralMultiplicity();
0552     internalId_.chgEMfrac_ = patjet->chargedEmEnergy() / jet->energy();
0553     internalId_.neuEMfrac_ = patjet->neutralEmEnergy() / jet->energy();
0554     internalId_.chgHadrfrac_ = patjet->chargedHadronEnergy() / jet->energy();
0555     internalId_.neuHadrfrac_ = patjet->neutralHadronEnergy() / jet->energy();
0556     if (usePuppi)
0557       internalId_.nNeutrals_ = multNeut;
0558   } else {
0559     internalId_.nCharged_ = pfjet->chargedMultiplicity();
0560     internalId_.nNeutrals_ = pfjet->neutralMultiplicity();
0561     internalId_.chgEMfrac_ = pfjet->chargedEmEnergy() / jet->energy();
0562     internalId_.neuEMfrac_ = pfjet->neutralEmEnergy() / jet->energy();
0563     internalId_.chgHadrfrac_ = pfjet->chargedHadronEnergy() / jet->energy();
0564     internalId_.neuHadrfrac_ = pfjet->neutralHadronEnergy() / jet->energy();
0565   }
0566   internalId_.nParticles_ = jet->nConstituents();
0567 
0568   ///////////////////////pull variable///////////////////////////////////
0569   float sumW2(0.0);
0570   float sum_deta(0.0), sum_dphi(0.0);
0571   float ave_deta(0.0), ave_dphi(0.0);
0572   for (size_t j = 0; j < jet->numberOfDaughters(); j++) {
0573     const auto& part = jet->daughterPtr(j);
0574     if (!(part.isAvailable() && part.isNonnull())) {
0575       continue;
0576     }
0577 
0578     float partPuppiWeight = 1.0;
0579     if (usePuppi) {
0580       const pat::PackedCandidate* partpack = dynamic_cast<const pat::PackedCandidate*>(part.get());
0581       if (partpack != nullptr)
0582         partPuppiWeight = partpack->puppiWeight();
0583     }
0584 
0585     float weight = (part->pt()) * partPuppiWeight;
0586     float weight2 = weight * weight;
0587     sumW2 += weight2;
0588     float deta = part->eta() - jet->eta();
0589     float dphi = reco::deltaPhi(*part, *jet);
0590     sum_deta += deta * weight2;
0591     sum_dphi += dphi * weight2;
0592     if (sumW2 > 0) {
0593       ave_deta = sum_deta / sumW2;
0594       ave_dphi = sum_dphi / sumW2;
0595     }
0596   }
0597 
0598   float ddetaR_sum(0.0), ddphiR_sum(0.0), pull_tmp(0.0);
0599   for (size_t i = 0; i < jet->numberOfDaughters(); i++) {
0600     const auto& part = jet->daughterPtr(i);
0601     if (!(part.isAvailable() && part.isNonnull())) {
0602       continue;
0603     }
0604 
0605     float partPuppiWeight = 1.0;
0606     if (usePuppi) {
0607       const pat::PackedCandidate* partpack = dynamic_cast<const pat::PackedCandidate*>(part.get());
0608       if (partpack != nullptr)
0609         partPuppiWeight = partpack->puppiWeight();
0610     }
0611 
0612     float weight = partPuppiWeight * (part->pt()) * partPuppiWeight * (part->pt());
0613     float deta = part->eta() - jet->eta();
0614     float dphi = reco::deltaPhi(*part, *jet);
0615     float ddeta, ddphi, ddR;
0616     ddeta = deta - ave_deta;
0617     ddphi = dphi - ave_dphi;
0618     ddR = sqrt(ddeta * ddeta + ddphi * ddphi);
0619     ddetaR_sum += ddR * ddeta * weight;
0620     ddphiR_sum += ddR * ddphi * weight;
0621   }
0622   if (sumW2 > 0) {
0623     float ddetaR_ave = ddetaR_sum / sumW2;
0624     float ddphiR_ave = ddphiR_sum / sumW2;
0625     pull_tmp = sqrt(ddetaR_ave * ddetaR_ave + ddphiR_ave * ddphiR_ave);
0626   }
0627   internalId_.pull_ = pull_tmp;
0628   ///////////////////////////////////////////////////////////////////////
0629 
0630   setPtEtaPhi(*lLead, internalId_.leadPt_, internalId_.leadEta_, internalId_.leadPhi_);
0631   setPtEtaPhi(*lSecond, internalId_.secondPt_, internalId_.secondEta_, internalId_.secondPhi_);
0632   setPtEtaPhi(*lLeadNeut, internalId_.leadNeutPt_, internalId_.leadNeutEta_, internalId_.leadNeutPhi_);
0633   setPtEtaPhi(*lLeadEm, internalId_.leadEmPt_, internalId_.leadEmEta_, internalId_.leadEmPhi_);
0634   setPtEtaPhi(*lLeadCh, internalId_.leadChPt_, internalId_.leadChEta_, internalId_.leadChPhi_);
0635 
0636   std::sort(frac.begin(), frac.end(), std::greater<float>());
0637   std::sort(fracCh.begin(), fracCh.end(), std::greater<float>());
0638   std::sort(fracEm.begin(), fracEm.end(), std::greater<float>());
0639   std::sort(fracNeut.begin(), fracNeut.end(), std::greater<float>());
0640   assign(frac, internalId_.leadFrac_, internalId_.secondFrac_, internalId_.thirdFrac_, internalId_.fourthFrac_);
0641   assign(
0642       fracCh, internalId_.leadChFrac_, internalId_.secondChFrac_, internalId_.thirdChFrac_, internalId_.fourthChFrac_);
0643   assign(
0644       fracEm, internalId_.leadEmFrac_, internalId_.secondEmFrac_, internalId_.thirdEmFrac_, internalId_.fourthEmFrac_);
0645   assign(fracNeut,
0646          internalId_.leadNeutFrac_,
0647          internalId_.secondNeutFrac_,
0648          internalId_.thirdNeutFrac_,
0649          internalId_.fourthNeutFrac_);
0650 
0651   covMatrix(0, 0) /= sumPt2;
0652   covMatrix(0, 1) /= sumPt2;
0653   covMatrix(1, 1) /= sumPt2;
0654   covMatrix(1, 0) = covMatrix(0, 1);
0655   internalId_.etaW_ = sqrt(covMatrix(0, 0));
0656   internalId_.phiW_ = sqrt(covMatrix(1, 1));
0657   internalId_.jetW_ = 0.5 * (internalId_.etaW_ + internalId_.phiW_);
0658   TVectorD eigVals(2);
0659   eigVals = TMatrixDSymEigen(covMatrix).GetEigenValues();
0660   internalId_.majW_ = sqrt(std::abs(eigVals(0)));
0661   internalId_.minW_ = sqrt(std::abs(eigVals(1)));
0662   if (internalId_.majW_ < internalId_.minW_) {
0663     std::swap(internalId_.majW_, internalId_.minW_);
0664   }
0665 
0666   internalId_.dRLeadCent_ = reco::deltaR(*jet, *lLead);
0667   if (lSecond != nullptr) {
0668     internalId_.dRLead2nd_ = reco::deltaR(*jet, *lSecond);
0669   }
0670   internalId_.dRMean_ /= jetPt;
0671   internalId_.dRMeanNeut_ /= jetPt;
0672   internalId_.dRMeanEm_ /= jetPt;
0673   internalId_.dRMeanCh_ /= jetPt;
0674   internalId_.dR2Mean_ /= sumPt2;
0675 
0676   for (size_t ic = 0; ic < ncones; ++ic) {
0677     *coneFracs[ic] /= jetPt;
0678     *coneEmFracs[ic] /= jetPt;
0679     *coneNeutFracs[ic] /= jetPt;
0680     *coneChFracs[ic] /= jetPt;
0681   }
0682   //http://jets.physics.harvard.edu/qvg/
0683   double ptMean = sumPt / internalId_.nParticles_;
0684   double ptRMS = 0;
0685   for (unsigned int i0 = 0; i0 < frac.size(); i0++) {
0686     ptRMS += (frac[i0] - ptMean) * (frac[i0] - ptMean);
0687   }
0688   ptRMS /= internalId_.nParticles_;
0689   ptRMS = sqrt(ptRMS);
0690 
0691   internalId_.ptMean_ = ptMean;
0692   internalId_.ptRMS_ = ptRMS / jetPt;
0693   internalId_.pt2A_ = sqrt(internalId_.ptD_ / internalId_.nParticles_) / jetPt;
0694   internalId_.ptD_ = sqrt(internalId_.ptD_) / sumPt;
0695   internalId_.ptDCh_ = sqrt(internalId_.ptDCh_) / sumPtCh;
0696   internalId_.ptDNe_ = sqrt(internalId_.ptDNe_) / sumPtNe;
0697   internalId_.sumPt_ = sumPt;
0698   internalId_.sumChPt_ = sumPtCh;
0699   internalId_.sumNePt_ = sumPtNe;
0700 
0701   internalId_.jetR_ = lLead->pt() / sumPt;
0702   internalId_.jetRchg_ = lLeadCh->pt() / sumPt;
0703   internalId_.dRMatch_ = dRmin;
0704 
0705   if (sumTkPt != 0.) {
0706     internalId_.beta_ /= sumTkPt;
0707     internalId_.betaStar_ /= sumTkPt;
0708     internalId_.betaClassic_ /= sumTkPt;
0709     internalId_.betaStarClassic_ /= sumTkPt;
0710   } else {
0711     assert(internalId_.beta_ == 0. && internalId_.betaStar_ == 0. && internalId_.betaClassic_ == 0. &&
0712            internalId_.betaStarClassic_ == 0.);
0713   }
0714 
0715   if (cache_->runMvas()) {
0716     runMva();
0717   }
0718 
0719   return PileupJetIdentifier(internalId_);
0720 }
0721 
0722 // ------------------------------------------------------------------------------------------
0723 std::string PileupJetIdAlgo::dumpVariables() const {
0724   std::stringstream out;
0725   for (variables_list_t::const_iterator it = variables_.begin(); it != variables_.end(); ++it) {
0726     out << std::setw(15) << it->first << std::setw(3) << "=" << std::setw(5) << *it->second.first << " ("
0727         << std::setw(5) << it->second.second << ")" << std::endl;
0728   }
0729   return out.str();
0730 }
0731 
0732 // ------------------------------------------------------------------------------------------
0733 void PileupJetIdAlgo::resetVariables() {
0734   internalId_.idFlag_ = 0;
0735   for (variables_list_t::iterator it = variables_.begin(); it != variables_.end(); ++it) {
0736     *it->second.first = it->second.second;
0737   }
0738 }
0739 
0740 // ------------------------------------------------------------------------------------------
0741 #define INIT_VARIABLE(NAME, TMVANAME, VAL) \
0742   internalId_.NAME##_ = VAL;               \
0743   variables_[#NAME] = std::make_pair(&internalId_.NAME##_, VAL);
0744 
0745 // ------------------------------------------------------------------------------------------
0746 void PileupJetIdAlgo::initVariables() {
0747   internalId_.idFlag_ = 0;
0748   INIT_VARIABLE(mva, "", -100.);
0749   //INIT_VARIABLE(jetPt      , "jspt_1", 0.);
0750   //INIT_VARIABLE(jetEta     , "jseta_1", large_val);
0751   INIT_VARIABLE(jetPt, "", 0.);
0752   INIT_VARIABLE(jetEta, "", large_val);
0753   INIT_VARIABLE(jetPhi, "jsphi_1", large_val);
0754   INIT_VARIABLE(jetM, "jm_1", 0.);
0755 
0756   INIT_VARIABLE(nCharged, "", 0.);
0757   INIT_VARIABLE(nNeutrals, "", 0.);
0758 
0759   INIT_VARIABLE(chgEMfrac, "", 0.);
0760   INIT_VARIABLE(neuEMfrac, "", 0.);
0761   INIT_VARIABLE(chgHadrfrac, "", 0.);
0762   INIT_VARIABLE(neuHadrfrac, "", 0.);
0763 
0764   INIT_VARIABLE(d0, "jd0_1", -1000.);
0765   INIT_VARIABLE(dZ, "jdZ_1", -1000.);
0766   //INIT_VARIABLE(nParticles , "npart_1"  , 0.);
0767   INIT_VARIABLE(nParticles, "", 0.);
0768 
0769   INIT_VARIABLE(leadPt, "lpt_1", 0.);
0770   INIT_VARIABLE(leadEta, "leta_1", large_val);
0771   INIT_VARIABLE(leadPhi, "lphi_1", large_val);
0772   INIT_VARIABLE(secondPt, "spt_1", 0.);
0773   INIT_VARIABLE(secondEta, "seta_1", large_val);
0774   INIT_VARIABLE(secondPhi, "sphi_1", large_val);
0775   INIT_VARIABLE(leadNeutPt, "lnept_1", 0.);
0776   INIT_VARIABLE(leadNeutEta, "lneeta_1", large_val);
0777   INIT_VARIABLE(leadNeutPhi, "lnephi_1", large_val);
0778   INIT_VARIABLE(leadEmPt, "lempt_1", 0.);
0779   INIT_VARIABLE(leadEmEta, "lemeta_1", large_val);
0780   INIT_VARIABLE(leadEmPhi, "lemphi_1", large_val);
0781   INIT_VARIABLE(leadChPt, "lchpt_1", 0.);
0782   INIT_VARIABLE(leadChEta, "lcheta_1", large_val);
0783   INIT_VARIABLE(leadChPhi, "lchphi_1", large_val);
0784   INIT_VARIABLE(leadFrac, "lLfr_1", 0.);
0785 
0786   INIT_VARIABLE(dRLeadCent, "drlc_1", 0.);
0787   INIT_VARIABLE(dRLead2nd, "drls_1", 0.);
0788   INIT_VARIABLE(dRMean, "drm_1", 0.);
0789   INIT_VARIABLE(dRMean, "", 0.);
0790   INIT_VARIABLE(pull, "", 0.);
0791   INIT_VARIABLE(dRMeanNeut, "drmne_1", 0.);
0792   INIT_VARIABLE(dRMeanEm, "drem_1", 0.);
0793   INIT_VARIABLE(dRMeanCh, "drch_1", 0.);
0794   INIT_VARIABLE(dR2Mean, "", 0.);
0795 
0796   INIT_VARIABLE(ptD, "", 0.);
0797   INIT_VARIABLE(ptMean, "", 0.);
0798   INIT_VARIABLE(ptRMS, "", 0.);
0799   INIT_VARIABLE(pt2A, "", 0.);
0800   INIT_VARIABLE(ptDCh, "", 0.);
0801   INIT_VARIABLE(ptDNe, "", 0.);
0802   INIT_VARIABLE(sumPt, "", 0.);
0803   INIT_VARIABLE(sumChPt, "", 0.);
0804   INIT_VARIABLE(sumNePt, "", 0.);
0805 
0806   INIT_VARIABLE(secondFrac, "", 0.);
0807   INIT_VARIABLE(thirdFrac, "", 0.);
0808   INIT_VARIABLE(fourthFrac, "", 0.);
0809 
0810   INIT_VARIABLE(leadChFrac, "", 0.);
0811   INIT_VARIABLE(secondChFrac, "", 0.);
0812   INIT_VARIABLE(thirdChFrac, "", 0.);
0813   INIT_VARIABLE(fourthChFrac, "", 0.);
0814 
0815   INIT_VARIABLE(leadNeutFrac, "", 0.);
0816   INIT_VARIABLE(secondNeutFrac, "", 0.);
0817   INIT_VARIABLE(thirdNeutFrac, "", 0.);
0818   INIT_VARIABLE(fourthNeutFrac, "", 0.);
0819 
0820   INIT_VARIABLE(leadEmFrac, "", 0.);
0821   INIT_VARIABLE(secondEmFrac, "", 0.);
0822   INIT_VARIABLE(thirdEmFrac, "", 0.);
0823   INIT_VARIABLE(fourthEmFrac, "", 0.);
0824 
0825   INIT_VARIABLE(jetW, "", 1.);
0826   INIT_VARIABLE(etaW, "", 1.);
0827   INIT_VARIABLE(phiW, "", 1.);
0828 
0829   INIT_VARIABLE(majW, "", 1.);
0830   INIT_VARIABLE(minW, "", 1.);
0831 
0832   INIT_VARIABLE(frac01, "", 0.);
0833   INIT_VARIABLE(frac02, "", 0.);
0834   INIT_VARIABLE(frac03, "", 0.);
0835   INIT_VARIABLE(frac04, "", 0.);
0836   INIT_VARIABLE(frac05, "", 0.);
0837   INIT_VARIABLE(frac06, "", 0.);
0838   INIT_VARIABLE(frac07, "", 0.);
0839 
0840   INIT_VARIABLE(chFrac01, "", 0.);
0841   INIT_VARIABLE(chFrac02, "", 0.);
0842   INIT_VARIABLE(chFrac03, "", 0.);
0843   INIT_VARIABLE(chFrac04, "", 0.);
0844   INIT_VARIABLE(chFrac05, "", 0.);
0845   INIT_VARIABLE(chFrac06, "", 0.);
0846   INIT_VARIABLE(chFrac07, "", 0.);
0847 
0848   INIT_VARIABLE(neutFrac01, "", 0.);
0849   INIT_VARIABLE(neutFrac02, "", 0.);
0850   INIT_VARIABLE(neutFrac03, "", 0.);
0851   INIT_VARIABLE(neutFrac04, "", 0.);
0852   INIT_VARIABLE(neutFrac05, "", 0.);
0853   INIT_VARIABLE(neutFrac06, "", 0.);
0854   INIT_VARIABLE(neutFrac07, "", 0.);
0855 
0856   INIT_VARIABLE(emFrac01, "", 0.);
0857   INIT_VARIABLE(emFrac02, "", 0.);
0858   INIT_VARIABLE(emFrac03, "", 0.);
0859   INIT_VARIABLE(emFrac04, "", 0.);
0860   INIT_VARIABLE(emFrac05, "", 0.);
0861   INIT_VARIABLE(emFrac06, "", 0.);
0862   INIT_VARIABLE(emFrac07, "", 0.);
0863 
0864   INIT_VARIABLE(beta, "", 0.);
0865   INIT_VARIABLE(betaStar, "", 0.);
0866   INIT_VARIABLE(betaClassic, "", 0.);
0867   INIT_VARIABLE(betaStarClassic, "", 0.);
0868 
0869   INIT_VARIABLE(nvtx, "", 0.);
0870   INIT_VARIABLE(rho, "", 0.);
0871   INIT_VARIABLE(nTrueInt, "", 0.);
0872 
0873   INIT_VARIABLE(jetR, "", 0.);
0874   INIT_VARIABLE(jetRchg, "", 0.);
0875   INIT_VARIABLE(dRMatch, "", 0.);
0876 }
0877 
0878 #undef INIT_VARIABLE