File indexing completed on 2025-03-28 03:57:23
0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/Frameworkfwd.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/Framework/interface/stream/EDProducer.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/Utilities/interface/StreamID.h"
0007 #include "RecoBTag/FeatureTools/interface/TrackInfoBuilder.h"
0008
0009 using namespace btagbtvdeep;
0010
0011 #include <string>
0012
0013 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0014 #include "DataFormats/Math/interface/deltaR.h"
0015 #include "DataFormats/NanoAOD/interface/FlatTable.h"
0016
0017
0018 #include "DataFormats/BTauReco/interface/JetTag.h"
0019 #include "DataFormats/BTauReco/interface/UnifiedParticleTransformerAK4TagInfo.h"
0020 #include "FWCore/Framework/interface/makeRefToBaseProdFrom.h"
0021
0022
0023 #include "DataFormats/PatCandidates/interface/Jet.h"
0024
0025 template <typename T>
0026 class JetTaggerTableProducer : public edm::stream::EDProducer<> {
0027 public:
0028 explicit JetTaggerTableProducer(const edm::ParameterSet &);
0029 ~JetTaggerTableProducer() override;
0030
0031 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0032
0033 private:
0034 void produce(edm::Event &, const edm::EventSetup &) override;
0035
0036 const std::string nameDeepJet_;
0037 const std::string idx_nameDeepJet_;
0038 const unsigned int n_cpf_ = 29;
0039 const unsigned int n_npf_ = 25;
0040 const unsigned int n_sv_ = 12;
0041 const unsigned int n_lt_ = 5;
0042
0043 edm::EDGetTokenT<edm::View<T>> jet_token_;
0044
0045 typedef std::vector<reco::UnifiedParticleTransformerAK4TagInfo> TagInfoCollection;
0046 const edm::EDGetTokenT<TagInfoCollection> tag_info_src_;
0047
0048 constexpr static bool usePhysForLightAndUndefined = false;
0049 };
0050
0051
0052
0053
0054 template <typename T>
0055 JetTaggerTableProducer<T>::JetTaggerTableProducer(const edm::ParameterSet &iConfig)
0056 : nameDeepJet_(iConfig.getParameter<std::string>("nameDeepJet")),
0057 idx_nameDeepJet_(iConfig.getParameter<std::string>("idx_nameDeepJet")),
0058 n_cpf_(iConfig.getParameter<unsigned int>("n_cpf")),
0059 n_npf_(iConfig.getParameter<unsigned int>("n_npf")),
0060 n_sv_(iConfig.getParameter<unsigned int>("n_sv")),
0061 n_lt_(iConfig.getParameter<unsigned int>("n_lt")),
0062 jet_token_(consumes<edm::View<T>>(iConfig.getParameter<edm::InputTag>("jets"))),
0063 tag_info_src_(consumes<TagInfoCollection>(iConfig.getParameter<edm::InputTag>("tagInfo_src"))) {
0064 produces<nanoaod::FlatTable>(nameDeepJet_);
0065 }
0066
0067 template <typename T>
0068 JetTaggerTableProducer<T>::~JetTaggerTableProducer() {}
0069
0070 template <typename T>
0071 void JetTaggerTableProducer<T>::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
0072
0073
0074
0075
0076 auto jets = iEvent.getHandle(jet_token_);
0077
0078 edm::Handle<TagInfoCollection> tag_infos;
0079 iEvent.getByToken(tag_info_src_, tag_infos);
0080
0081 unsigned nJets = jets->size();
0082
0083 std::vector<int> jet_N_CPFCands(nJets);
0084 std::vector<int> jet_N_NPFCands(nJets);
0085 std::vector<int> jet_N_SVs(nJets);
0086 std::vector<int> jet_N_LTs(nJets);
0087
0088
0089 std::vector<std::vector<float>> Cpfcan_BtagPf_trackEtaRel_nCpf(n_cpf_, std::vector<float>(nJets));
0090 std::vector<std::vector<float>> Cpfcan_BtagPf_trackPtRel_nCpf(n_cpf_, std::vector<float>(nJets));
0091 std::vector<std::vector<float>> Cpfcan_BtagPf_trackPPar_nCpf(n_cpf_, std::vector<float>(nJets));
0092 std::vector<std::vector<float>> Cpfcan_BtagPf_trackDeltaR_nCpf(n_cpf_, std::vector<float>(nJets));
0093 std::vector<std::vector<float>> Cpfcan_BtagPf_trackPParRatio_nCpf(n_cpf_, std::vector<float>(nJets));
0094 std::vector<std::vector<float>> Cpfcan_BtagPf_trackSip2dVal_nCpf(n_cpf_, std::vector<float>(nJets));
0095 std::vector<std::vector<float>> Cpfcan_BtagPf_trackSip2dSig_nCpf(n_cpf_, std::vector<float>(nJets));
0096 std::vector<std::vector<float>> Cpfcan_BtagPf_trackSip3dVal_nCpf(n_cpf_, std::vector<float>(nJets));
0097 std::vector<std::vector<float>> Cpfcan_BtagPf_trackSip3dSig_nCpf(n_cpf_, std::vector<float>(nJets));
0098 std::vector<std::vector<float>> Cpfcan_BtagPf_trackJetDistVal_nCpf(n_cpf_, std::vector<float>(nJets));
0099 std::vector<std::vector<float>> Cpfcan_ptrel_nCpf(n_cpf_, std::vector<float>(nJets));
0100 std::vector<std::vector<float>> Cpfcan_drminsv_nCpf(n_cpf_, std::vector<float>(nJets));
0101 std::vector<std::vector<int>> Cpfcan_VTX_ass_nCpf(n_cpf_, std::vector<int>(nJets));
0102 std::vector<std::vector<float>> Cpfcan_puppiw_nCpf(n_cpf_, std::vector<float>(nJets));
0103 std::vector<std::vector<float>> Cpfcan_chi2_nCpf(n_cpf_, std::vector<float>(nJets));
0104 std::vector<std::vector<int>> Cpfcan_quality_nCpf(n_cpf_, std::vector<int>(nJets));
0105 std::vector<std::vector<float>> Cpfcan_charge_nCpf(n_cpf_, std::vector<float>(nJets));
0106 std::vector<std::vector<float>> Cpfcan_dz_nCpf(n_cpf_, std::vector<float>(nJets));
0107 std::vector<std::vector<float>> Cpfcan_btagPf_trackDecayLen_nCpf(n_cpf_, std::vector<float>(nJets));
0108 std::vector<std::vector<float>> Cpfcan_HadFrac_nCpf(n_cpf_, std::vector<float>(nJets));
0109 std::vector<std::vector<float>> Cpfcan_CaloFrac_nCpf(n_cpf_, std::vector<float>(nJets));
0110 std::vector<std::vector<float>> Cpfcan_pdgID_nCpf(n_cpf_, std::vector<float>(nJets));
0111 std::vector<std::vector<float>> Cpfcan_lostInnerHits_nCpf(n_cpf_, std::vector<float>(nJets));
0112 std::vector<std::vector<float>> Cpfcan_numberOfPixelHits_nCpf(n_cpf_, std::vector<float>(nJets));
0113 std::vector<std::vector<float>> Cpfcan_numberOfStripHits_nCpf(n_cpf_, std::vector<float>(nJets));
0114 std::vector<std::vector<float>> Cpfcan_px_nCpf(n_cpf_, std::vector<float>(nJets));
0115 std::vector<std::vector<float>> Cpfcan_py_nCpf(n_cpf_, std::vector<float>(nJets));
0116 std::vector<std::vector<float>> Cpfcan_pz_nCpf(n_cpf_, std::vector<float>(nJets));
0117 std::vector<std::vector<float>> Cpfcan_e_nCpf(n_cpf_, std::vector<float>(nJets));
0118
0119
0120 std::vector<std::vector<float>> Npfcan_ptrel_nNpf(n_npf_, std::vector<float>(nJets));
0121 std::vector<std::vector<float>> Npfcan_etarel_nNpf(n_npf_, std::vector<float>(nJets));
0122 std::vector<std::vector<float>> Npfcan_phirel_nNpf(n_npf_, std::vector<float>(nJets));
0123 std::vector<std::vector<float>> Npfcan_deltaR_nNpf(n_npf_, std::vector<float>(nJets));
0124 std::vector<std::vector<int>> Npfcan_isGamma_nNpf(n_npf_, std::vector<int>(nJets));
0125 std::vector<std::vector<float>> Npfcan_HadFrac_nNpf(n_npf_, std::vector<float>(nJets));
0126 std::vector<std::vector<float>> Npfcan_drminsv_nNpf(n_npf_, std::vector<float>(nJets));
0127 std::vector<std::vector<float>> Npfcan_puppiw_nNpf(n_npf_, std::vector<float>(nJets));
0128 std::vector<std::vector<float>> Npfcan_px_nNpf(n_npf_, std::vector<float>(nJets));
0129 std::vector<std::vector<float>> Npfcan_py_nNpf(n_npf_, std::vector<float>(nJets));
0130 std::vector<std::vector<float>> Npfcan_pz_nNpf(n_npf_, std::vector<float>(nJets));
0131 std::vector<std::vector<float>> Npfcan_e_nNpf(n_npf_, std::vector<float>(nJets));
0132
0133
0134 std::vector<std::vector<float>> sv_pt_nSV(n_sv_, std::vector<float>(nJets));
0135 std::vector<std::vector<float>> sv_deltaR_nSV(n_sv_, std::vector<float>(nJets));
0136 std::vector<std::vector<float>> sv_mass_nSV(n_sv_, std::vector<float>(nJets));
0137 std::vector<std::vector<float>> sv_etarel_nSV(n_sv_, std::vector<float>(nJets));
0138 std::vector<std::vector<float>> sv_phirel_nSV(n_sv_, std::vector<float>(nJets));
0139 std::vector<std::vector<int>> sv_ntracks_nSV(n_sv_, std::vector<int>(nJets));
0140 std::vector<std::vector<float>> sv_chi2_nSV(n_sv_, std::vector<float>(nJets));
0141 std::vector<std::vector<float>> sv_normchi2_nSV(n_sv_, std::vector<float>(nJets));
0142 std::vector<std::vector<float>> sv_dxy_nSV(n_sv_, std::vector<float>(nJets));
0143 std::vector<std::vector<float>> sv_dxysig_nSV(n_sv_, std::vector<float>(nJets));
0144 std::vector<std::vector<float>> sv_d3d_nSV(n_sv_, std::vector<float>(nJets));
0145 std::vector<std::vector<float>> sv_d3dsig_nSV(n_sv_, std::vector<float>(nJets));
0146 std::vector<std::vector<float>> sv_costhetasvpv_nSV(n_sv_, std::vector<float>(nJets));
0147 std::vector<std::vector<float>> sv_enratio_nSV(n_sv_, std::vector<float>(nJets));
0148 #ifdef JTTP_NEED_SV_PE
0149 std::vector<std::vector<float>> sv_px_nSV(n_sv_, std::vector<float>(nJets));
0150 std::vector<std::vector<float>> sv_py_nSV(n_sv_, std::vector<float>(nJets));
0151 std::vector<std::vector<float>> sv_pz_nSV(n_sv_, std::vector<float>(nJets));
0152 std::vector<std::vector<float>> sv_e_nSV(n_sv_, std::vector<float>(nJets));
0153 #else
0154 std::vector<std::vector<float>> sv_eta_nSV(n_sv_, std::vector<float>(nJets));
0155 std::vector<std::vector<float>> sv_phi_nSV(n_sv_, std::vector<float>(nJets));
0156 #endif
0157
0158
0159 std::vector<std::vector<float>> lt_btagPf_trackEtaRel_nLT(n_lt_, std::vector<float>(nJets));
0160 std::vector<std::vector<float>> lt_btagPf_trackPtRel_nLT(n_lt_, std::vector<float>(nJets));
0161 std::vector<std::vector<float>> lt_btagPf_trackPPar_nLT(n_lt_, std::vector<float>(nJets));
0162 std::vector<std::vector<float>> lt_btagPf_trackDeltaR_nLT(n_lt_, std::vector<float>(nJets));
0163 std::vector<std::vector<float>> lt_btagPf_trackPParRatio_nLT(n_lt_, std::vector<float>(nJets));
0164 std::vector<std::vector<float>> lt_btagPf_trackSip2dVal_nLT(n_lt_, std::vector<float>(nJets));
0165 std::vector<std::vector<float>> lt_btagPf_trackSip2dSig_nLT(n_lt_, std::vector<float>(nJets));
0166 std::vector<std::vector<float>> lt_btagPf_trackSip3dVal_nLT(n_lt_, std::vector<float>(nJets));
0167 std::vector<std::vector<float>> lt_btagPf_trackSip3dSig_nLT(n_lt_, std::vector<float>(nJets));
0168 std::vector<std::vector<float>> lt_btagPf_trackJetDistVal_nLT(n_lt_, std::vector<float>(nJets));
0169 std::vector<std::vector<float>> lt_drminsv_nLT(n_lt_, std::vector<float>(nJets));
0170 std::vector<std::vector<float>> lt_charge_nLT(n_lt_, std::vector<float>(nJets));
0171 std::vector<std::vector<float>> lt_puppiw_nLT(n_lt_, std::vector<float>(nJets));
0172 std::vector<std::vector<float>> lt_chi2_nLT(n_lt_, std::vector<float>(nJets));
0173 std::vector<std::vector<float>> lt_quality_nLT(n_lt_, std::vector<float>(nJets));
0174 std::vector<std::vector<float>> lt_lostInnerHits_nLT(n_lt_, std::vector<float>(nJets));
0175 std::vector<std::vector<float>> lt_numberOfPixelHits_nLT(n_lt_, std::vector<float>(nJets));
0176 std::vector<std::vector<float>> lt_numberOfStripHits_nLT(n_lt_, std::vector<float>(nJets));
0177 std::vector<std::vector<float>> lt_pt_nLT(n_lt_, std::vector<float>(nJets));
0178 std::vector<std::vector<float>> lt_eta_nLT(n_lt_, std::vector<float>(nJets));
0179 std::vector<std::vector<float>> lt_phi_nLT(n_lt_, std::vector<float>(nJets));
0180 std::vector<std::vector<float>> lt_e_nLT(n_lt_, std::vector<float>(nJets));
0181
0182 if (!tag_infos->empty()) {
0183 for (unsigned i_jet = 0; i_jet < nJets; ++i_jet) {
0184
0185
0186 const auto &taginfo = (*tag_infos)[i_jet];
0187 const auto &features = taginfo.features();
0188
0189
0190
0191 jet_N_CPFCands[i_jet] = features.c_pf_features.size();
0192 jet_N_NPFCands[i_jet] = features.n_pf_features.size();
0193 jet_N_SVs[i_jet] = features.sv_features.size();
0194 jet_N_LTs[i_jet] = features.lt_features.size();
0195
0196 std::vector<const btagbtvdeep::ChargedCandidateFeatures *> ranked_c_pf_features;
0197 ranked_c_pf_features.reserve(features.c_pf_features.size());
0198 for (auto &c_pf : features.c_pf_features)
0199 ranked_c_pf_features.push_back(&c_pf);
0200
0201 std::vector<const btagbtvdeep::NeutralCandidateFeatures *> ranked_n_pf_features;
0202 ranked_n_pf_features.reserve(features.n_pf_features.size());
0203 for (auto &n_pf : features.n_pf_features)
0204 ranked_n_pf_features.push_back(&n_pf);
0205
0206 std::vector<const btagbtvdeep::SecondaryVertexFeatures *> ranked_sv_features;
0207 ranked_sv_features.reserve(features.sv_features.size());
0208 for (auto &sv : features.sv_features)
0209 ranked_sv_features.push_back(&sv);
0210
0211 std::vector<const btagbtvdeep::LostTracksFeatures *> ranked_lt_features;
0212 ranked_lt_features.reserve(features.lt_features.size());
0213 for (auto < : features.lt_features)
0214 ranked_lt_features.push_back(<);
0215
0216 auto max_c_pf_n = std::min(features.c_pf_features.size(), (std::size_t)n_cpf_);
0217 auto max_n_pf_n = std::min(features.n_pf_features.size(), (std::size_t)n_npf_);
0218 auto max_sv_n = std::min(features.sv_features.size(), (std::size_t)n_sv_);
0219 auto max_lt_n = std::min(features.lt_features.size(), (std::size_t)n_lt_);
0220
0221 auto c_pf_cmp = [](const btagbtvdeep::ChargedCandidateFeatures *a,
0222 const btagbtvdeep::ChargedCandidateFeatures *b) { return a->pt > b->pt; };
0223
0224
0225 auto sv_cmp = [](const btagbtvdeep::SecondaryVertexFeatures *a, const btagbtvdeep::SecondaryVertexFeatures *b) {
0226 return a->pt > b->pt;
0227 };
0228
0229
0230
0231 auto c_pf_cmp_ip = [](const btagbtvdeep::ChargedCandidateFeatures *a,
0232 const btagbtvdeep::ChargedCandidateFeatures *b) {
0233 return fabs(a->btagPf_trackSip3dVal) > fabs(b->btagPf_trackSip3dVal);
0234 };
0235 auto sv_cmp_ip = [](const btagbtvdeep::SecondaryVertexFeatures *a,
0236 const btagbtvdeep::SecondaryVertexFeatures *b) { return fabs(a->d3d) > fabs(b->d3d); };
0237
0238
0239 if (n_cpf_ == 2) {
0240
0241
0242 if (!ranked_c_pf_features.empty()) {
0243 auto highest_pT = *std::min_element(ranked_c_pf_features.begin(), ranked_c_pf_features.end(), c_pf_cmp);
0244 auto highest_IP = *std::min_element(ranked_c_pf_features.begin(), ranked_c_pf_features.end(), c_pf_cmp_ip);
0245 ranked_c_pf_features = {highest_pT, highest_IP};
0246 }
0247 } else {
0248
0249
0250
0251
0252 }
0253
0254
0255 {
0256
0257
0258
0259
0260 }
0261
0262
0263 if (n_sv_ == 2) {
0264
0265
0266 if (ranked_sv_features.size() >= 2) {
0267 auto highest_pT = *std::min_element(ranked_sv_features.begin(), ranked_sv_features.end(), sv_cmp);
0268 auto highest_IP = *std::min_element(ranked_sv_features.begin(), ranked_sv_features.end(), sv_cmp_ip);
0269 if (highest_IP == highest_pT) {
0270 std::nth_element(
0271 ranked_sv_features.begin(), next(ranked_sv_features.begin()), ranked_sv_features.end(), sv_cmp_ip);
0272 for (size_t isv = 0; isv < 2; ++isv) {
0273 highest_IP = ranked_sv_features[isv];
0274 if (highest_IP != highest_pT)
0275 break;
0276 }
0277 }
0278 ranked_sv_features.clear();
0279 ranked_sv_features = {highest_pT, highest_IP};
0280 }
0281 } else {
0282
0283
0284
0285
0286 }
0287
0288
0289 {
0290
0291
0292
0293
0294 }
0295
0296
0297 for (std::size_t c_pf_n = 0; c_pf_n < max_c_pf_n; c_pf_n++) {
0298 const auto &c_pf_features = *ranked_c_pf_features.at(c_pf_n);
0299 Cpfcan_BtagPf_trackEtaRel_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackEtaRel;
0300 Cpfcan_BtagPf_trackPtRel_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackPtRel;
0301 Cpfcan_BtagPf_trackPPar_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackPPar;
0302 Cpfcan_BtagPf_trackDeltaR_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackDeltaR;
0303 Cpfcan_BtagPf_trackPParRatio_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackPParRatio;
0304 Cpfcan_BtagPf_trackSip2dVal_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackSip2dVal;
0305 Cpfcan_BtagPf_trackSip2dSig_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackSip2dSig;
0306 Cpfcan_BtagPf_trackSip3dVal_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackSip3dVal;
0307 Cpfcan_BtagPf_trackSip3dSig_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackSip3dSig;
0308 Cpfcan_BtagPf_trackJetDistVal_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackJetDistVal;
0309 Cpfcan_ptrel_nCpf[c_pf_n][i_jet] = c_pf_features.ptrel;
0310 Cpfcan_drminsv_nCpf[c_pf_n][i_jet] = c_pf_features.drminsv;
0311 Cpfcan_VTX_ass_nCpf[c_pf_n][i_jet] = c_pf_features.vtx_ass;
0312 Cpfcan_puppiw_nCpf[c_pf_n][i_jet] = c_pf_features.puppiw;
0313 Cpfcan_chi2_nCpf[c_pf_n][i_jet] = c_pf_features.chi2;
0314 Cpfcan_quality_nCpf[c_pf_n][i_jet] = c_pf_features.quality;
0315 Cpfcan_charge_nCpf[c_pf_n][i_jet] = c_pf_features.charge;
0316 Cpfcan_dz_nCpf[c_pf_n][i_jet] = c_pf_features.dz;
0317 Cpfcan_btagPf_trackDecayLen_nCpf[c_pf_n][i_jet] = c_pf_features.btagPf_trackDecayLen;
0318 Cpfcan_HadFrac_nCpf[c_pf_n][i_jet] = c_pf_features.HadFrac;
0319 Cpfcan_CaloFrac_nCpf[c_pf_n][i_jet] = c_pf_features.CaloFrac;
0320 Cpfcan_pdgID_nCpf[c_pf_n][i_jet] = c_pf_features.pdgID;
0321 Cpfcan_lostInnerHits_nCpf[c_pf_n][i_jet] = c_pf_features.lostInnerHits;
0322 Cpfcan_numberOfPixelHits_nCpf[c_pf_n][i_jet] = c_pf_features.numberOfPixelHits;
0323 Cpfcan_numberOfStripHits_nCpf[c_pf_n][i_jet] = c_pf_features.numberOfStripHits;
0324 Cpfcan_px_nCpf[c_pf_n][i_jet] = c_pf_features.px;
0325 Cpfcan_py_nCpf[c_pf_n][i_jet] = c_pf_features.py;
0326 Cpfcan_pz_nCpf[c_pf_n][i_jet] = c_pf_features.pz;
0327 Cpfcan_e_nCpf[c_pf_n][i_jet] = c_pf_features.e;
0328 }
0329
0330
0331 for (std::size_t n_pf_n = 0; n_pf_n < max_n_pf_n; n_pf_n++) {
0332 const auto &n_pf_features = *ranked_n_pf_features.at(n_pf_n);
0333 Npfcan_ptrel_nNpf[n_pf_n][i_jet] = n_pf_features.ptrel;
0334 Npfcan_etarel_nNpf[n_pf_n][i_jet] = n_pf_features.etarel;
0335 Npfcan_phirel_nNpf[n_pf_n][i_jet] = n_pf_features.phirel;
0336 Npfcan_deltaR_nNpf[n_pf_n][i_jet] = n_pf_features.deltaR;
0337 Npfcan_isGamma_nNpf[n_pf_n][i_jet] = n_pf_features.isGamma;
0338 Npfcan_HadFrac_nNpf[n_pf_n][i_jet] = n_pf_features.hadFrac;
0339 Npfcan_drminsv_nNpf[n_pf_n][i_jet] = n_pf_features.drminsv;
0340 Npfcan_puppiw_nNpf[n_pf_n][i_jet] = n_pf_features.puppiw;
0341 Npfcan_px_nNpf[n_pf_n][i_jet] = n_pf_features.px;
0342 Npfcan_py_nNpf[n_pf_n][i_jet] = n_pf_features.py;
0343 Npfcan_pz_nNpf[n_pf_n][i_jet] = n_pf_features.pz;
0344 Npfcan_e_nNpf[n_pf_n][i_jet] = n_pf_features.e;
0345 }
0346
0347
0348 for (std::size_t sv_n = 0; sv_n < max_sv_n; sv_n++) {
0349 const auto &sv_features = *ranked_sv_features.at(sv_n);
0350 sv_pt_nSV[sv_n][i_jet] = sv_features.pt;
0351 sv_deltaR_nSV[sv_n][i_jet] = sv_features.deltaR;
0352 sv_mass_nSV[sv_n][i_jet] = sv_features.mass;
0353 sv_etarel_nSV[sv_n][i_jet] = sv_features.etarel;
0354 sv_phirel_nSV[sv_n][i_jet] = sv_features.phirel;
0355 sv_ntracks_nSV[sv_n][i_jet] = sv_features.ntracks;
0356 sv_chi2_nSV[sv_n][i_jet] = sv_features.chi2;
0357 sv_normchi2_nSV[sv_n][i_jet] = sv_features.normchi2;
0358 sv_dxy_nSV[sv_n][i_jet] = sv_features.dxy;
0359 sv_dxysig_nSV[sv_n][i_jet] = sv_features.dxysig;
0360 sv_d3d_nSV[sv_n][i_jet] = sv_features.d3d;
0361 sv_d3dsig_nSV[sv_n][i_jet] = sv_features.d3dsig;
0362 sv_costhetasvpv_nSV[sv_n][i_jet] = sv_features.costhetasvpv;
0363 sv_enratio_nSV[sv_n][i_jet] = sv_features.enratio;
0364 #ifdef JTTP_NEED_SV_PE
0365 sv_px_nSV[sv_n][i_jet] = sv_features.px;
0366 sv_py_nSV[sv_n][i_jet] = sv_features.py;
0367 sv_pz_nSV[sv_n][i_jet] = sv_features.pz;
0368 sv_e_nSV[sv_n][i_jet] = sv_features.e;
0369 #else
0370 sv_eta_nSV[sv_n][i_jet] = sv_features.eta;
0371 sv_phi_nSV[sv_n][i_jet] = sv_features.phi;
0372 #endif
0373 }
0374
0375
0376 for (std::size_t lt_n = 0; lt_n < max_lt_n; lt_n++) {
0377 const auto <_features = *ranked_lt_features.at(lt_n);
0378 lt_btagPf_trackEtaRel_nLT[lt_n][i_jet] = lt_features.btagPf_trackEtaRel;
0379 lt_btagPf_trackPtRel_nLT[lt_n][i_jet] = lt_features.btagPf_trackPtRel;
0380 lt_btagPf_trackPPar_nLT[lt_n][i_jet] = lt_features.btagPf_trackPPar;
0381 lt_btagPf_trackDeltaR_nLT[lt_n][i_jet] = lt_features.btagPf_trackDeltaR;
0382 lt_btagPf_trackPParRatio_nLT[lt_n][i_jet] = lt_features.btagPf_trackPParRatio;
0383 lt_btagPf_trackSip2dVal_nLT[lt_n][i_jet] = lt_features.btagPf_trackSip2dVal;
0384 lt_btagPf_trackSip2dSig_nLT[lt_n][i_jet] = lt_features.btagPf_trackSip2dSig;
0385 lt_btagPf_trackSip3dVal_nLT[lt_n][i_jet] = lt_features.btagPf_trackSip3dVal;
0386 lt_btagPf_trackSip3dSig_nLT[lt_n][i_jet] = lt_features.btagPf_trackSip3dSig;
0387 lt_btagPf_trackJetDistVal_nLT[lt_n][i_jet] = lt_features.btagPf_trackJetDistVal;
0388 lt_drminsv_nLT[lt_n][i_jet] = lt_features.drminsv;
0389 lt_charge_nLT[lt_n][i_jet] = lt_features.charge;
0390 lt_puppiw_nLT[lt_n][i_jet] = lt_features.puppiw;
0391 lt_chi2_nLT[lt_n][i_jet] = lt_features.chi2;
0392 lt_quality_nLT[lt_n][i_jet] = lt_features.quality;
0393 lt_lostInnerHits_nLT[lt_n][i_jet] = lt_features.lostInnerHits;
0394 lt_numberOfPixelHits_nLT[lt_n][i_jet] = lt_features.numberOfPixelHits;
0395 lt_numberOfStripHits_nLT[lt_n][i_jet] = lt_features.numberOfStripHits;
0396 lt_pt_nLT[lt_n][i_jet] = lt_features.pt;
0397 lt_eta_nLT[lt_n][i_jet] = lt_features.eta;
0398 lt_phi_nLT[lt_n][i_jet] = lt_features.phi;
0399 lt_e_nLT[lt_n][i_jet] = lt_features.e;
0400 }
0401 }
0402 }
0403
0404
0405 auto djTable = std::make_unique<nanoaod::FlatTable>(jet_N_CPFCands.size(), nameDeepJet_, false, true);
0406
0407 djTable->addColumn<int>("DeepJet_nCpfcand", jet_N_CPFCands, "Number of charged PF candidates in the jet");
0408 djTable->addColumn<int>("DeepJet_nNpfcand", jet_N_NPFCands, "Number of neutral PF candidates in the jet");
0409 djTable->addColumn<int>("DeepJet_nsv", jet_N_SVs, "Number of secondary vertices in the jet");
0410 djTable->addColumn<int>("DeepJet_nlt", jet_N_LTs, "Number of lost tracks in the jet");
0411
0412
0413 for (unsigned int p = 0; p < n_cpf_; p++) {
0414 auto s = std::to_string(p);
0415
0416 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackDeltaR_" + s,
0417 Cpfcan_BtagPf_trackDeltaR_nCpf[p],
0418 "track pseudoangular distance from the jet axis for the " + s + ". cpf",
0419 10);
0420 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackEtaRel_" + s,
0421 Cpfcan_BtagPf_trackEtaRel_nCpf[p],
0422 "track pseudorapidity, relative to the jet axis for the " + s + ". cpf",
0423 10);
0424 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackJetDistVal_" + s,
0425 Cpfcan_BtagPf_trackJetDistVal_nCpf[p],
0426 "minimum track approach distance to jet axis for the " + s + ". cpf",
0427 10);
0428 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackPPar_" + s,
0429 Cpfcan_BtagPf_trackPPar_nCpf[p],
0430 "dot product of the jet and track momentum for the " + s + ". cpf",
0431 10);
0432 djTable->addColumn<float>(
0433 "DeepJet_Cpfcan_BtagPf_trackPParRatio_" + s,
0434 Cpfcan_BtagPf_trackPParRatio_nCpf[p],
0435 "dot product of the jet and track momentum divided by the magnitude of the jet momentum for the " + s + ". cpf",
0436 10);
0437 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackPtRel_" + s,
0438 Cpfcan_BtagPf_trackPtRel_nCpf[p],
0439 "track transverse momentum, relative to the jet axis for the " + s + ". cpf",
0440 10);
0441 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackSip2dSig_" + s,
0442 Cpfcan_BtagPf_trackSip2dSig_nCpf[p],
0443 "track 2D signed impact parameter significance for the " + s + ". cpf",
0444 10);
0445 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackSip3dSig_" + s,
0446 Cpfcan_BtagPf_trackSip3dSig_nCpf[p],
0447 "track 3D signed impact parameter significance for the " + s + ". cpf",
0448 10);
0449 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackSip2dVal_" + s,
0450 Cpfcan_BtagPf_trackSip2dVal_nCpf[p],
0451 "track 2D signed impact parameter for the " + s + ". cpf",
0452 10);
0453 djTable->addColumn<float>("DeepJet_Cpfcan_BtagPf_trackSip3dVal_" + s,
0454 Cpfcan_BtagPf_trackSip3dVal_nCpf[p],
0455 "track 3D signed impact parameter for the " + s + ". cpf",
0456 10);
0457 djTable->addColumn<float>("DeepJet_Cpfcan_ptrel_" + s,
0458 Cpfcan_ptrel_nCpf[p],
0459 "fraction of the jet momentum carried by the track for the " + s + ". cpf",
0460 10);
0461 djTable->addColumn<float>("DeepJet_Cpfcan_drminsv_" + s,
0462 Cpfcan_drminsv_nCpf[p],
0463 "track pseudoangular distance from the closest secondary vertex of the " + s + ". cpf",
0464 10);
0465 djTable->addColumn<int>(
0466 "DeepJet_Cpfcan_VTX_ass_" + s,
0467 Cpfcan_VTX_ass_nCpf[p],
0468 "integer flag that indicates whether the track was used in the primary vertex fit for the " + s + ". cpf",
0469 10);
0470 djTable->addColumn<float>("DeepJet_Cpfcan_puppiw_" + s,
0471 Cpfcan_puppiw_nCpf[p],
0472 "charged candidate PUPPI weight of the " + s + ". cpf",
0473 10);
0474 djTable->addColumn<float>(
0475 "DeepJet_Cpfcan_chi2_" + s, Cpfcan_chi2_nCpf[p], "chi2 of the charged track fit for the " + s + ". cpf", 10);
0476 djTable->addColumn<int>(
0477 "DeepJet_Cpfcan_quality_" + s,
0478 Cpfcan_quality_nCpf[p],
0479 "integer flag which indicates the quality of the fitted track, based on number of detector hits used for the "
0480 "reconstruction as well as the overall chi2 of the charged track fit for the " +
0481 s + ". cpf",
0482 10);
0483 djTable->addColumn<float>("DeepJet_Cpfcan_charge_" + s, Cpfcan_charge_nCpf[p], "", 10);
0484 djTable->addColumn<float>("DeepJet_Cpfcan_dz_" + s, Cpfcan_dz_nCpf[p], "", 10);
0485 djTable->addColumn<float>("DeepJet_Cpfcan_btagPf_trackDecayLen_" + s, Cpfcan_btagPf_trackDecayLen_nCpf[p], "", 10);
0486 djTable->addColumn<float>("DeepJet_Cpfcan_HadFrac_" + s, Cpfcan_HadFrac_nCpf[p], "", 10);
0487 djTable->addColumn<float>("DeepJet_Cpfcan_CaloFrac_" + s, Cpfcan_CaloFrac_nCpf[p], "", 10);
0488 djTable->addColumn<float>("DeepJet_Cpfcan_pdgID_" + s, Cpfcan_pdgID_nCpf[p], "", 10);
0489 djTable->addColumn<float>("DeepJet_Cpfcan_lostInnerHits_" + s, Cpfcan_lostInnerHits_nCpf[p], "", 10);
0490 djTable->addColumn<float>("DeepJet_Cpfcan_numberOfPixelHits_" + s, Cpfcan_numberOfPixelHits_nCpf[p], "", 10);
0491 djTable->addColumn<float>("DeepJet_Cpfcan_numberOfStripHits_" + s, Cpfcan_numberOfStripHits_nCpf[p], "", 10);
0492 djTable->addColumn<float>("DeepJet_Cpfcan_px_" + s, Cpfcan_px_nCpf[p], "", 10);
0493 djTable->addColumn<float>("DeepJet_Cpfcan_py_" + s, Cpfcan_py_nCpf[p], "", 10);
0494 djTable->addColumn<float>("DeepJet_Cpfcan_pz_" + s, Cpfcan_pz_nCpf[p], "", 10);
0495 djTable->addColumn<float>("DeepJet_Cpfcan_e_" + s, Cpfcan_e_nCpf[p], "", 10);
0496 }
0497
0498
0499 for (unsigned int p = 0; p < n_npf_; p++) {
0500 auto s = std::to_string(p);
0501
0502 djTable->addColumn<float>("DeepJet_Npfcan_ptrel_" + s,
0503 Npfcan_ptrel_nNpf[p],
0504 "fraction of the jet momentum carried by the neutral candidate for the " + s + ". npf",
0505 10);
0506 djTable->addColumn<float>("DeepJetExtra_Npfcan_etarel_" + s,
0507 Npfcan_etarel_nNpf[p],
0508 "pseudorapidity relative to parent jet for the " + s + ". npf",
0509 10);
0510 djTable->addColumn<float>(
0511 "DeepJetExtra_Npfcan_phirel_" + s, Npfcan_phirel_nNpf[p], "DeltaPhi(npf, jet) for the " + s + ". npf", 10);
0512 djTable->addColumn<float>(
0513 "DeepJet_Npfcan_deltaR_" + s,
0514 Npfcan_deltaR_nNpf[p],
0515 "pseudoangular distance between the neutral candidate and the jet axis for the " + s + ". npf",
0516 10);
0517 djTable->addColumn<int>("DeepJet_Npfcan_isGamma_" + s,
0518 Npfcan_isGamma_nNpf[p],
0519 "integer flag indicating whether the neutral candidate is a photon for the " + s + ". npf",
0520 10);
0521 djTable->addColumn<float>(
0522 "DeepJet_Npfcan_HadFrac_" + s,
0523 Npfcan_HadFrac_nNpf[p],
0524 "fraction of the neutral candidate energy deposited in the hadronic calorimeter for the " + s + ". npf",
0525 10);
0526 djTable->addColumn<float>(
0527 "DeepJet_Npfcan_drminsv_" + s,
0528 Npfcan_drminsv_nNpf[p],
0529 "pseudoangular distance between the neutral candidate and the closest secondary vertex for the " + s + ". npf",
0530 10);
0531 djTable->addColumn<float>("DeepJet_Npfcan_puppiw_" + s,
0532 Npfcan_puppiw_nNpf[p],
0533 "neutral candidate PUPPI weight for the " + s + ". npf",
0534 10);
0535 djTable->addColumn<float>("DeepJet_Npfcan_px_" + s, Npfcan_px_nNpf[p], "", 10);
0536 djTable->addColumn<float>("DeepJet_Npfcan_py_" + s, Npfcan_py_nNpf[p], "", 10);
0537 djTable->addColumn<float>("DeepJet_Npfcan_pz_" + s, Npfcan_pz_nNpf[p], "", 10);
0538 djTable->addColumn<float>("DeepJet_Npfcan_e_" + s, Npfcan_e_nNpf[p], "", 10);
0539 }
0540
0541
0542 for (unsigned int p = 0; p < n_sv_; p++) {
0543 auto s = std::to_string(p);
0544
0545 djTable->addColumn<float>("DeepJet_sv_pt_" + s, sv_pt_nSV[p], "SV pt of the " + s + ". SV", 10);
0546 djTable->addColumn<float>("DeepJet_sv_deltaR_" + s,
0547 sv_deltaR_nSV[p],
0548 "pseudoangular distance between jet axis and the " + s + ". SV direction",
0549 10);
0550 djTable->addColumn<float>("DeepJet_sv_mass_" + s, sv_mass_nSV[p], "SV mass of the " + s + ". SV", 10);
0551 djTable->addColumn<float>("DeepJetExtra_sv_etarel_" + s,
0552 sv_etarel_nSV[p],
0553 "pseudorapidity relative to parent jet for the " + s + ". SV",
0554 10);
0555 djTable->addColumn<float>(
0556 "DeepJetExtra_sv_phirel_" + s, sv_phirel_nSV[p], "DeltaPhi(sv, jet) for the " + s + ". SV", 10);
0557 djTable->addColumn<float>(
0558 "DeepJet_sv_ntracks_" + s, sv_ntracks_nSV[p], "Number of tracks asociated to the " + s + ". SV", 10);
0559 djTable->addColumn<float>("DeepJet_sv_chi2_" + s, sv_chi2_nSV[p], "chi2 of the " + s + ". SV", 10);
0560 djTable->addColumn<float>("DeepJet_sv_normchi2_" + s, sv_normchi2_nSV[p], "chi2/dof of the " + s + ". SV", 10);
0561 djTable->addColumn<float>(
0562 "DeepJet_sv_dxy_" + s, sv_dxy_nSV[p], "2D impact parameter (flight distance) value of the " + s + ". SV", 10);
0563 djTable->addColumn<float>("DeepJet_sv_dxysig_" + s,
0564 sv_dxysig_nSV[p],
0565 "2D impact parameter (flight distance) significance of the " + s + ". SV",
0566 10);
0567 djTable->addColumn<float>(
0568 "DeepJet_sv_d3d_" + s, sv_d3d_nSV[p], "3D impact parameter (flight distance) value of the " + s + ". SV", 10);
0569 djTable->addColumn<float>("DeepJet_sv_d3dsig_" + s,
0570 sv_d3dsig_nSV[p],
0571 "3D impact parameter (flight distance) significance of the " + s + ". SV",
0572 10);
0573 djTable->addColumn<float>("DeepJet_sv_costhetasvpv_" + s,
0574 sv_costhetasvpv_nSV[p],
0575 "cosine of the angle between the " + s +
0576 ". SV flight direction and the direction of the " + s + ". SV momentum",
0577 10);
0578 djTable->addColumn<float>(
0579 "DeepJet_sv_enratio_" + s, sv_enratio_nSV[p], "ratio of the " + s + ". SV energy ratio to the jet energy", 10);
0580 #ifdef JTTP_NEED_SV_PE
0581 djTable->addColumn<float>("DeepJet_sv_px_" + s, sv_px_nSV[p], "", 10);
0582 djTable->addColumn<float>("DeepJet_sv_py_" + s, sv_py_nSV[p], "", 10);
0583 djTable->addColumn<float>("DeepJet_sv_pz_" + s, sv_pz_nSV[p], "", 10);
0584 djTable->addColumn<float>("DeepJet_sv_e_" + s, sv_e_nSV[p], "", 10);
0585 #else
0586 djTable->addColumn<float>("DeepJet_sv_eta_" + s, sv_eta_nSV[p], "", 10);
0587 djTable->addColumn<float>("DeepJet_sv_phi_" + s, sv_phi_nSV[p], "", 10);
0588 #endif
0589 }
0590
0591
0592 for (unsigned int p = 0; p < n_lt_; p++) {
0593 auto s = std::to_string(p);
0594
0595 djTable->addColumn<float>("DeepJet_lt_btagPf_trackEtaRel_" + s, lt_btagPf_trackEtaRel_nLT[p], "", 10);
0596 djTable->addColumn<float>("DeepJet_lt_btagPf_trackPtRel_" + s, lt_btagPf_trackPtRel_nLT[p], "", 10);
0597 djTable->addColumn<float>("DeepJet_lt_btagPf_trackPPar_" + s, lt_btagPf_trackPPar_nLT[p], "", 10);
0598 djTable->addColumn<float>("DeepJet_lt_btagPf_trackDeltaR_" + s, lt_btagPf_trackDeltaR_nLT[p], "", 10);
0599 djTable->addColumn<float>("DeepJet_lt_btagPf_trackPParRatio_" + s, lt_btagPf_trackPParRatio_nLT[p], "", 10);
0600 djTable->addColumn<float>("DeepJet_lt_btagPf_trackSip2dVal_" + s, lt_btagPf_trackSip2dVal_nLT[p], "", 10);
0601 djTable->addColumn<float>("DeepJet_lt_btagPf_trackSip2dSig_" + s, lt_btagPf_trackSip2dSig_nLT[p], "", 10);
0602 djTable->addColumn<float>("DeepJet_lt_btagPf_trackSip3dVal_" + s, lt_btagPf_trackSip3dVal_nLT[p], "", 10);
0603 djTable->addColumn<float>("DeepJet_lt_btagPf_trackSip3dSig_" + s, lt_btagPf_trackSip3dSig_nLT[p], "", 10);
0604 djTable->addColumn<float>("DeepJet_lt_btagPf_trackJetDistVal_" + s, lt_btagPf_trackJetDistVal_nLT[p], "", 10);
0605 djTable->addColumn<float>("DeepJet_lt_drminsv_" + s, lt_drminsv_nLT[p], "", 10);
0606 djTable->addColumn<float>("DeepJet_lt_charge_" + s, lt_charge_nLT[p], "", 10);
0607 djTable->addColumn<float>("DeepJet_lt_puppiw_" + s, lt_puppiw_nLT[p], "", 10);
0608 djTable->addColumn<float>("DeepJet_lt_chi2_" + s, lt_chi2_nLT[p], "", 10);
0609 djTable->addColumn<float>("DeepJet_lt_quality_" + s, lt_quality_nLT[p], "", 10);
0610 djTable->addColumn<float>("DeepJet_lt_lostInnerHits_" + s, lt_lostInnerHits_nLT[p], "", 10);
0611 djTable->addColumn<float>("DeepJet_lt_numberOfPixelHits_" + s, lt_numberOfPixelHits_nLT[p], "", 10);
0612 djTable->addColumn<float>("DeepJet_lt_numberOfStripHits_" + s, lt_numberOfStripHits_nLT[p], "", 10);
0613 djTable->addColumn<float>("DeepJet_lt_pt_" + s, lt_pt_nLT[p], "", 10);
0614 djTable->addColumn<float>("DeepJet_lt_eta_" + s, lt_eta_nLT[p], "", 10);
0615 djTable->addColumn<float>("DeepJet_lt_phi_" + s, lt_phi_nLT[p], "", 10);
0616 djTable->addColumn<float>("DeepJet_lt_e_" + s, lt_e_nLT[p], "", 10);
0617 }
0618
0619 iEvent.put(std::move(djTable), nameDeepJet_);
0620 }
0621
0622 template <typename T>
0623 void JetTaggerTableProducer<T>::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0624 edm::ParameterSetDescription desc;
0625 desc.add<std::string>("nameDeepJet", "Jet");
0626 desc.add<std::string>("idx_nameDeepJet", "djIdx");
0627
0628 desc.add<unsigned int>("n_cpf", 2);
0629 desc.add<unsigned int>("n_npf", 2);
0630 desc.add<unsigned int>("n_sv", 2);
0631 desc.add<unsigned int>("n_lt", 2);
0632 desc.add<edm::InputTag>("jets", edm::InputTag("slimmedJetsPuppi"));
0633 desc.add<edm::InputTag>("tagInfo_src", edm::InputTag("pfUnifiedParticleTransformerAK4TagInfosPuppiWithDeepInfo"));
0634 descriptions.addWithDefaultLabel(desc);
0635 }
0636
0637 typedef JetTaggerTableProducer<pat::Jet> PatJetTaggerTableProducer;
0638
0639 DEFINE_FWK_MODULE(PatJetTaggerTableProducer);