Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:26

0001 
0002 #include "FWCore/Framework/interface/Frameworkfwd.h"
0003 #include "FWCore/Framework/interface/stream/EDProducer.h"
0004 
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/Utilities/interface/StreamID.h"
0010 #include "FWCore/Utilities/interface/ESGetToken.h"
0011 
0012 #include "DataFormats/PatCandidates/interface/Jet.h"
0013 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0014 
0015 #include "DataFormats/BTauReco/interface/ShallowTagInfo.h"
0016 
0017 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0018 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0019 
0020 #include "DataFormats/BTauReco/interface/DeepFlavourTagInfo.h"
0021 #include "DataFormats/BTauReco/interface/DeepFlavourFeatures.h"
0022 
0023 #include "RecoBTag/FeatureTools/interface/JetConverter.h"
0024 #include "RecoBTag/FeatureTools/interface/ShallowTagInfoConverter.h"
0025 #include "RecoBTag/FeatureTools/interface/SecondaryVertexConverter.h"
0026 #include "RecoBTag/FeatureTools/interface/NeutralCandidateConverter.h"
0027 #include "RecoBTag/FeatureTools/interface/ChargedCandidateConverter.h"
0028 
0029 #include "RecoBTag/FeatureTools/interface/TrackInfoBuilder.h"
0030 #include "RecoBTag/FeatureTools/interface/sorting_modules.h"
0031 
0032 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0033 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0034 
0035 #include "RecoBTag/FeatureTools/interface/deep_helpers.h"
0036 
0037 #include "FWCore/ParameterSet/interface/Registry.h"
0038 #include "FWCore/Common/interface/Provenance.h"
0039 #include "DataFormats/Provenance/interface/ProductProvenance.h"
0040 
0041 #include "DataFormats/BTauReco/interface/SeedingTrackFeatures.h"
0042 #include "DataFormats/BTauReco/interface/TrackPairFeatures.h"
0043 #include "RecoBTag/FeatureTools/interface/TrackPairInfoBuilder.h"
0044 #include "RecoBTag/FeatureTools/interface/SeedingTrackInfoBuilder.h"
0045 #include "RecoBTag/FeatureTools/interface/SeedingTracksConverter.h"
0046 
0047 #include "FWCore/Framework/interface/ESHandle.h"
0048 #include "FWCore/Framework/interface/EventSetup.h"
0049 #include "RecoBTag/TrackProbability/interface/HistogramProbabilityEstimator.h"
0050 class HistogramProbabilityEstimator;
0051 #include <memory>
0052 
0053 #include <typeinfo>
0054 #include "CondFormats/BTauObjects/interface/TrackProbabilityCalibration.h"
0055 #include "CondFormats/DataRecord/interface/BTagTrackProbability2DRcd.h"
0056 #include "CondFormats/DataRecord/interface/BTagTrackProbability3DRcd.h"
0057 #include "FWCore/Framework/interface/EventSetupRecord.h"
0058 #include "FWCore/Framework/interface/EventSetupRecordImplementation.h"
0059 #include "FWCore/Framework/interface/EventSetupRecordKey.h"
0060 #include "DataFormats/Common/interface/AssociationMap.h"
0061 
0062 class DeepFlavourTagInfoProducer : public edm::stream::EDProducer<> {
0063 public:
0064   explicit DeepFlavourTagInfoProducer(const edm::ParameterSet&);
0065   ~DeepFlavourTagInfoProducer() override;
0066 
0067   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0068 
0069 private:
0070   typedef std::vector<reco::DeepFlavourTagInfo> DeepFlavourTagInfoCollection;
0071   typedef reco::VertexCompositePtrCandidateCollection SVCollection;
0072   typedef reco::VertexCollection VertexCollection;
0073   typedef edm::View<reco::ShallowTagInfo> ShallowTagInfoCollection;
0074   typedef edm::AssociationMap<edm::OneToOne<reco::JetView, reco::JetView>> JetMatchMap;
0075 
0076   void beginStream(edm::StreamID) override {}
0077   void produce(edm::Event&, const edm::EventSetup&) override;
0078   void endStream() override {}
0079 
0080   const double jet_radius_;
0081   const double min_candidate_pt_;
0082   const bool flip_;
0083 
0084   edm::EDGetTokenT<edm::View<reco::Jet>> jet_token_;
0085   edm::EDGetTokenT<JetMatchMap> unsubjet_map_token_;
0086   edm::EDGetTokenT<VertexCollection> vtx_token_;
0087   edm::EDGetTokenT<SVCollection> sv_token_;
0088   edm::EDGetTokenT<ShallowTagInfoCollection> shallow_tag_info_token_;
0089   edm::EDGetTokenT<edm::ValueMap<float>> puppi_value_map_token_;
0090   edm::EDGetTokenT<edm::ValueMap<int>> pvasq_value_map_token_;
0091   edm::EDGetTokenT<edm::Association<VertexCollection>> pvas_token_;
0092   edm::EDGetTokenT<edm::View<reco::Candidate>> candidateToken_;
0093   edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> track_builder_token_;
0094   edm::ESGetToken<TrackProbabilityCalibration, BTagTrackProbability2DRcd> calib2d_token_;
0095   edm::ESGetToken<TrackProbabilityCalibration, BTagTrackProbability3DRcd> calib3d_token_;
0096 
0097   bool use_puppi_value_map_;
0098   bool use_pvasq_value_map_;
0099   bool use_unsubjet_map_;
0100 
0101   bool fallback_puppi_weight_;
0102   bool fallback_vertex_association_;
0103 
0104   bool run_deepVertex_;
0105 
0106   bool is_weighted_jet_;
0107 
0108   //TrackProbability
0109   void checkEventSetup(const edm::EventSetup& iSetup);
0110   std::unique_ptr<HistogramProbabilityEstimator> probabilityEstimator_;
0111   bool compute_probabilities_;
0112   unsigned long long calibrationCacheId2D_;
0113   unsigned long long calibrationCacheId3D_;
0114 
0115   const double min_jet_pt_;
0116   const double max_jet_eta_;
0117 };
0118 
0119 DeepFlavourTagInfoProducer::DeepFlavourTagInfoProducer(const edm::ParameterSet& iConfig)
0120     : jet_radius_(iConfig.getParameter<double>("jet_radius")),
0121       min_candidate_pt_(iConfig.getParameter<double>("min_candidate_pt")),
0122       flip_(iConfig.getParameter<bool>("flip")),
0123       jet_token_(consumes<edm::View<reco::Jet>>(iConfig.getParameter<edm::InputTag>("jets"))),
0124       vtx_token_(consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
0125       sv_token_(consumes<SVCollection>(iConfig.getParameter<edm::InputTag>("secondary_vertices"))),
0126       shallow_tag_info_token_(
0127           consumes<ShallowTagInfoCollection>(iConfig.getParameter<edm::InputTag>("shallow_tag_infos"))),
0128       candidateToken_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("candidates"))),
0129       track_builder_token_(
0130           esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"))),
0131       use_puppi_value_map_(false),
0132       use_pvasq_value_map_(false),
0133       use_unsubjet_map_(false),
0134       fallback_puppi_weight_(iConfig.getParameter<bool>("fallback_puppi_weight")),
0135       fallback_vertex_association_(iConfig.getParameter<bool>("fallback_vertex_association")),
0136       run_deepVertex_(iConfig.getParameter<bool>("run_deepVertex")),
0137       is_weighted_jet_(iConfig.getParameter<bool>("is_weighted_jet")),
0138       compute_probabilities_(iConfig.getParameter<bool>("compute_probabilities")),
0139       min_jet_pt_(iConfig.getParameter<double>("min_jet_pt")),
0140       max_jet_eta_(iConfig.getParameter<double>("max_jet_eta")) {
0141   produces<DeepFlavourTagInfoCollection>();
0142 
0143   const auto& puppi_value_map_tag = iConfig.getParameter<edm::InputTag>("puppi_value_map");
0144   if (!puppi_value_map_tag.label().empty()) {
0145     puppi_value_map_token_ = consumes<edm::ValueMap<float>>(puppi_value_map_tag);
0146     use_puppi_value_map_ = true;
0147   } else if (is_weighted_jet_) {
0148     throw edm::Exception(edm::errors::Configuration,
0149                          "puppi_value_map is not set but jet is weighted. Must set puppi_value_map.");
0150   }
0151 
0152   const auto& pvas_tag = iConfig.getParameter<edm::InputTag>("vertex_associator");
0153   if (!pvas_tag.label().empty()) {
0154     pvasq_value_map_token_ = consumes<edm::ValueMap<int>>(pvas_tag);
0155     pvas_token_ = consumes<edm::Association<VertexCollection>>(pvas_tag);
0156     use_pvasq_value_map_ = true;
0157   }
0158   if (compute_probabilities_) {
0159     calib2d_token_ = esConsumes<TrackProbabilityCalibration, BTagTrackProbability2DRcd>();
0160     calib3d_token_ = esConsumes<TrackProbabilityCalibration, BTagTrackProbability3DRcd>();
0161   }
0162 
0163   const auto& unsubjet_map_tag = iConfig.getParameter<edm::InputTag>("unsubjet_map");
0164   if (!unsubjet_map_tag.label().empty()) {
0165     unsubjet_map_token_ = consumes<JetMatchMap>(unsubjet_map_tag);
0166     use_unsubjet_map_ = true;
0167   }
0168 }
0169 
0170 DeepFlavourTagInfoProducer::~DeepFlavourTagInfoProducer() {}
0171 
0172 void DeepFlavourTagInfoProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0173   // pfDeepFlavourTagInfos
0174   edm::ParameterSetDescription desc;
0175   desc.add<edm::InputTag>("shallow_tag_infos", edm::InputTag("pfDeepCSVTagInfos"));
0176   desc.add<double>("jet_radius", 0.4);
0177   desc.add<double>("min_candidate_pt", 0.95);
0178   desc.add<bool>("flip", false);
0179   desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
0180   desc.add<edm::InputTag>("puppi_value_map", edm::InputTag("puppi"));
0181   desc.add<edm::InputTag>("secondary_vertices", edm::InputTag("inclusiveCandidateSecondaryVertices"));
0182   desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
0183   desc.add<edm::InputTag>("unsubjet_map", {});
0184   desc.add<edm::InputTag>("candidates", edm::InputTag("packedPFCandidates"));
0185   desc.add<edm::InputTag>("vertex_associator", edm::InputTag("primaryVertexAssociation", "original"));
0186   desc.add<bool>("fallback_puppi_weight", false);
0187   desc.add<bool>("fallback_vertex_association", false);
0188   desc.add<bool>("run_deepVertex", false);
0189   desc.add<bool>("is_weighted_jet", false);
0190   desc.add<bool>("compute_probabilities", false);
0191   desc.add<double>("min_jet_pt", 15.0);
0192   desc.add<double>("max_jet_eta", 2.5);
0193   descriptions.add("pfDeepFlavourTagInfos", desc);
0194 }
0195 
0196 void DeepFlavourTagInfoProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0197   auto output_tag_infos = std::make_unique<DeepFlavourTagInfoCollection>();
0198   if (compute_probabilities_)
0199     checkEventSetup(iSetup);
0200 
0201   edm::Handle<edm::View<reco::Jet>> jets;
0202   iEvent.getByToken(jet_token_, jets);
0203 
0204   edm::Handle<JetMatchMap> unsubjet_map;
0205   if (use_unsubjet_map_)
0206     iEvent.getByToken(unsubjet_map_token_, unsubjet_map);
0207 
0208   edm::Handle<VertexCollection> vtxs;
0209   iEvent.getByToken(vtx_token_, vtxs);
0210   if (vtxs->empty()) {
0211     // produce empty TagInfos in case no primary vertex
0212     iEvent.put(std::move(output_tag_infos));
0213     return;  // exit event
0214   }
0215   // reference to primary vertex
0216   const auto& pv = vtxs->at(0);
0217 
0218   edm::Handle<edm::View<reco::Candidate>> tracks;
0219   iEvent.getByToken(candidateToken_, tracks);
0220 
0221   edm::Handle<SVCollection> svs;
0222   iEvent.getByToken(sv_token_, svs);
0223 
0224   edm::Handle<ShallowTagInfoCollection> shallow_tag_infos;
0225   iEvent.getByToken(shallow_tag_info_token_, shallow_tag_infos);
0226   double negative_cut = 0;  //used only with flip_
0227   if (flip_) {              //FIXME: Check if can do even less often than once per event
0228     const edm::Provenance* prov = shallow_tag_infos.provenance();
0229     const edm::ParameterSet& psetFromProvenance = edm::parameterSet(prov->stable(), iEvent.processHistory());
0230     negative_cut = ((psetFromProvenance.getParameter<edm::ParameterSet>("computer"))
0231                         .getParameter<edm::ParameterSet>("trackSelection"))
0232                        .getParameter<double>("sip3dSigMax");
0233   }
0234 
0235   edm::Handle<edm::ValueMap<float>> puppi_value_map;
0236   if (use_puppi_value_map_) {
0237     iEvent.getByToken(puppi_value_map_token_, puppi_value_map);
0238   }
0239 
0240   edm::Handle<edm::ValueMap<int>> pvasq_value_map;
0241   edm::Handle<edm::Association<VertexCollection>> pvas;
0242   if (use_pvasq_value_map_) {
0243     iEvent.getByToken(pvasq_value_map_token_, pvasq_value_map);
0244     iEvent.getByToken(pvas_token_, pvas);
0245   }
0246 
0247   edm::ESHandle<TransientTrackBuilder> track_builder = iSetup.getHandle(track_builder_token_);
0248 
0249   std::vector<reco::TransientTrack> selectedTracks;
0250   std::vector<float> masses;
0251 
0252   if (run_deepVertex_)  //make a collection of selected transient tracks for deepvertex outside of the jet loop
0253   {
0254     for (typename edm::View<reco::Candidate>::const_iterator track = tracks->begin(); track != tracks->end(); ++track) {
0255       unsigned int k = track - tracks->begin();
0256       if (track->bestTrack() != nullptr && track->pt() > 0.5) {
0257         if (std::fabs(track->vz() - pv.z()) < 0.5) {
0258           selectedTracks.push_back(track_builder->build(tracks->ptrAt(k)));
0259           masses.push_back(track->mass());
0260         }
0261       }
0262     }
0263   }
0264 
0265   for (std::size_t jet_n = 0; jet_n < jets->size(); jet_n++) {
0266     // create data containing structure
0267     btagbtvdeep::DeepFlavourFeatures features;
0268 
0269     // reco jet reference (use as much as possible)
0270     const auto& jet = jets->at(jet_n);
0271     // dynamical casting to pointers, null if not possible
0272     const auto* pf_jet = dynamic_cast<const reco::PFJet*>(&jet);
0273     const auto* pat_jet = dynamic_cast<const pat::Jet*>(&jet);
0274     edm::RefToBase<reco::Jet> jet_ref(jets, jet_n);
0275     const auto& unsubJet =
0276         (use_unsubjet_map_ && (*unsubjet_map)[jet_ref].isNonnull()) ? *(*unsubjet_map)[jet_ref] : jet;
0277     // TagInfoCollection not in an associative container so search for matchs
0278     const edm::View<reco::ShallowTagInfo>& taginfos = *shallow_tag_infos;
0279     edm::Ptr<reco::ShallowTagInfo> match;
0280     // Try first by 'same index'
0281     if ((jet_n < taginfos.size()) && (taginfos[jet_n].jet() == jet_ref)) {
0282       match = taginfos.ptrAt(jet_n);
0283     } else {
0284       // otherwise fail back to a simple search
0285       for (auto itTI = taginfos.begin(), edTI = taginfos.end(); itTI != edTI; ++itTI) {
0286         if (itTI->jet() == jet_ref) {
0287           match = taginfos.ptrAt(itTI - taginfos.begin());
0288           break;
0289         }
0290       }
0291     }
0292     reco::ShallowTagInfo tag_info;
0293     if (match.isNonnull()) {
0294       tag_info = *match;
0295     }  // will be default values otherwise
0296 
0297     // fill basic jet features
0298     btagbtvdeep::JetConverter::jetToFeatures(jet, features.jet_features);
0299 
0300     // fill number of pv
0301     features.npv = vtxs->size();
0302     math::XYZVector jet_dir = jet.momentum().Unit();
0303     GlobalVector jet_ref_track_dir(jet.px(), jet.py(), jet.pz());
0304 
0305     // fill features from ShallowTagInfo
0306     const auto& tag_info_vars = tag_info.taggingVariables();
0307     btagbtvdeep::bTagToFeatures(tag_info_vars, features.tag_info_features);
0308 
0309     // copy which will be sorted
0310     auto svs_sorted = *svs;
0311     // sort by dxy
0312     std::sort(svs_sorted.begin(), svs_sorted.end(), [&pv](const auto& sva, const auto& svb) {
0313       return btagbtvdeep::sv_vertex_comparator(sva, svb, pv);
0314     });
0315     // fill features from secondary vertices
0316     for (const auto& sv : svs_sorted) {
0317       if (reco::deltaR2(sv, jet_dir) > (jet_radius_ * jet_radius_))
0318         continue;
0319       else {
0320         features.sv_features.emplace_back();
0321         // in C++17 could just get from emplace_back output
0322         auto& sv_features = features.sv_features.back();
0323         btagbtvdeep::svToFeatures(sv, pv, jet, sv_features, flip_);
0324       }
0325     }
0326 
0327     // stuff required for dealing with pf candidates
0328 
0329     std::vector<btagbtvdeep::SortingClass<size_t>> c_sorted, n_sorted;
0330 
0331     // to cache the TrackInfo
0332     std::map<unsigned int, btagbtvdeep::TrackInfoBuilder> trackinfos;
0333 
0334     // unsorted reference to sv
0335     const auto& svs_unsorted = *svs;
0336     // fill collection, from DeepTNtuples plus some styling
0337     for (unsigned int i = 0; i < unsubJet.numberOfDaughters(); i++) {
0338       auto cand = unsubJet.daughter(i);
0339       if (cand) {
0340         // candidates under 950MeV (configurable) are not considered
0341         // might change if we use also white-listing
0342         if (cand->pt() < min_candidate_pt_)
0343           continue;
0344         if (cand->charge() != 0) {
0345           auto& trackinfo = trackinfos.emplace(i, track_builder).first->second;
0346           trackinfo.buildTrackInfo(cand, jet_dir, jet_ref_track_dir, pv);
0347           c_sorted.emplace_back(
0348               i, trackinfo.getTrackSip2dSig(), -btagbtvdeep::mindrsvpfcand(svs_unsorted, cand), cand->pt() / jet.pt());
0349         } else {
0350           n_sorted.emplace_back(i, -1, -btagbtvdeep::mindrsvpfcand(svs_unsorted, cand), cand->pt() / jet.pt());
0351         }
0352       }
0353     }
0354 
0355     // sort collections (open the black-box if you please)
0356     std::sort(c_sorted.begin(), c_sorted.end(), btagbtvdeep::SortingClass<std::size_t>::compareByABCInv);
0357     std::sort(n_sorted.begin(), n_sorted.end(), btagbtvdeep::SortingClass<std::size_t>::compareByABCInv);
0358 
0359     std::vector<size_t> c_sortedindices, n_sortedindices;
0360 
0361     // this puts 0 everywhere and the right position in ind
0362     c_sortedindices = btagbtvdeep::invertSortingVector(c_sorted);
0363     n_sortedindices = btagbtvdeep::invertSortingVector(n_sorted);
0364 
0365     // set right size to vectors
0366     features.c_pf_features.clear();
0367     features.c_pf_features.resize(c_sorted.size());
0368     features.n_pf_features.clear();
0369     features.n_pf_features.resize(n_sorted.size());
0370 
0371     for (unsigned int i = 0; i < unsubJet.numberOfDaughters(); i++) {
0372       // get pointer and check that is correct
0373       auto cand = dynamic_cast<const reco::Candidate*>(unsubJet.daughter(i));
0374       if (!cand)
0375         continue;
0376       // candidates under 950MeV are not considered
0377       // might change if we use also white-listing
0378       if (cand->pt() < 0.95)
0379         continue;
0380 
0381       auto packed_cand = dynamic_cast<const pat::PackedCandidate*>(cand);
0382       auto reco_cand = dynamic_cast<const reco::PFCandidate*>(cand);
0383 
0384       // need some edm::Ptr or edm::Ref if reco candidates
0385       reco::PFCandidatePtr reco_ptr;
0386       if (pf_jet) {
0387         reco_ptr = pf_jet->getPFConstituent(i);
0388       } else if (pat_jet && reco_cand) {
0389         reco_ptr = pat_jet->getPFConstituent(i);
0390       }
0391 
0392       reco::CandidatePtr cand_ptr;
0393       if (pat_jet) {
0394         cand_ptr = pat_jet->sourceCandidatePtr(i);
0395       }
0396 
0397       //
0398       // Access puppi weight from ValueMap.
0399       //
0400       float puppiw = 1.0;  // Set to fallback value
0401 
0402       if (reco_cand) {
0403         if (use_puppi_value_map_)
0404           puppiw = (*puppi_value_map)[reco_ptr];
0405         else if (!fallback_puppi_weight_) {
0406           throw edm::Exception(edm::errors::InvalidReference, "PUPPI value map missing")
0407               << "use fallback_puppi_weight option to use " << puppiw << " for reco_cand as default";
0408         }
0409       } else if (packed_cand) {
0410         if (use_puppi_value_map_)
0411           puppiw = (*puppi_value_map)[cand_ptr];
0412         else if (!fallback_puppi_weight_) {
0413           throw edm::Exception(edm::errors::InvalidReference, "PUPPI value map missing")
0414               << "use fallback_puppi_weight option to use " << puppiw << " for packed_cand as default";
0415         }
0416       } else {
0417         throw edm::Exception(edm::errors::InvalidReference)
0418             << "Cannot convert to either reco::PFCandidate or pat::PackedCandidate";
0419       }
0420 
0421       float drminpfcandsv = btagbtvdeep::mindrsvpfcand(svs_unsorted, cand);
0422 
0423       if (cand->charge() != 0) {
0424         // is charged candidate
0425         auto entry = c_sortedindices.at(i);
0426         // get cached track info
0427         auto& trackinfo = trackinfos.at(i);
0428         if (flip_ && (trackinfo.getTrackSip3dSig() > negative_cut)) {
0429           continue;
0430         }
0431         // get_ref to vector element
0432         auto& c_pf_features = features.c_pf_features.at(entry);
0433         // fill feature structure
0434         if (packed_cand) {
0435           btagbtvdeep::packedCandidateToFeatures(packed_cand,
0436                                                  jet,
0437                                                  trackinfo,
0438                                                  is_weighted_jet_,
0439                                                  drminpfcandsv,
0440                                                  static_cast<float>(jet_radius_),
0441                                                  puppiw,
0442                                                  c_pf_features,
0443                                                  flip_);
0444         } else if (reco_cand) {
0445           // get vertex association quality
0446           int pv_ass_quality = 0;  // fallback value
0447           if (use_pvasq_value_map_) {
0448             pv_ass_quality = (*pvasq_value_map)[reco_ptr];
0449           } else if (!fallback_vertex_association_) {
0450             throw edm::Exception(edm::errors::InvalidReference, "vertex association missing")
0451                 << "use fallback_vertex_association option to use" << pv_ass_quality
0452                 << "as default quality and closest dz PV as criteria";
0453           }
0454           // getting the PV as PackedCandidatesProducer
0455           // but using not the slimmed but original vertices
0456           auto ctrack = reco_cand->bestTrack();
0457           int pvi = -1;
0458           float dist = 1e99;
0459           for (size_t ii = 0; ii < vtxs->size(); ii++) {
0460             float dz = (ctrack) ? std::abs(ctrack->dz(((*vtxs)[ii]).position())) : 0;
0461             if (dz < dist) {
0462               pvi = ii;
0463               dist = dz;
0464             }
0465           }
0466           auto PV = reco::VertexRef(vtxs, pvi);
0467           if (use_pvasq_value_map_) {
0468             const reco::VertexRef& PV_orig = (*pvas)[reco_ptr];
0469             if (PV_orig.isNonnull())
0470               PV = reco::VertexRef(vtxs, PV_orig.key());
0471           }
0472           btagbtvdeep::recoCandidateToFeatures(reco_cand,
0473                                                jet,
0474                                                trackinfo,
0475                                                is_weighted_jet_,
0476                                                drminpfcandsv,
0477                                                static_cast<float>(jet_radius_),
0478                                                puppiw,
0479                                                pv_ass_quality,
0480                                                PV,
0481                                                c_pf_features,
0482                                                flip_);
0483         }
0484       } else {
0485         // is neutral candidate
0486         auto entry = n_sortedindices.at(i);
0487         // get_ref to vector element
0488         auto& n_pf_features = features.n_pf_features.at(entry);
0489         // fill feature structure
0490         if (packed_cand) {
0491           btagbtvdeep::packedCandidateToFeatures(
0492               packed_cand, jet, is_weighted_jet_, drminpfcandsv, static_cast<float>(jet_radius_), puppiw, n_pf_features);
0493         } else if (reco_cand) {
0494           btagbtvdeep::recoCandidateToFeatures(
0495               reco_cand, jet, is_weighted_jet_, drminpfcandsv, static_cast<float>(jet_radius_), puppiw, n_pf_features);
0496         }
0497       }
0498     }
0499 
0500     if (run_deepVertex_) {
0501       if (jet.pt() > min_jet_pt_ && std::fabs(jet.eta()) < max_jet_eta_)  //jet thresholds
0502         btagbtvdeep::seedingTracksToFeatures(selectedTracks,
0503                                              masses,
0504                                              jet,
0505                                              pv,
0506                                              probabilityEstimator_.get(),
0507                                              compute_probabilities_,
0508                                              features.seed_features);
0509     }
0510 
0511     output_tag_infos->emplace_back(features, jet_ref);
0512   }
0513 
0514   iEvent.put(std::move(output_tag_infos));
0515 }
0516 
0517 void DeepFlavourTagInfoProducer::checkEventSetup(const edm::EventSetup& iSetup) {
0518   using namespace edm;
0519   using namespace edm::eventsetup;
0520 
0521   const EventSetupRecord& re2D = iSetup.get<BTagTrackProbability2DRcd>();
0522   const EventSetupRecord& re3D = iSetup.get<BTagTrackProbability3DRcd>();
0523   unsigned long long cacheId2D = re2D.cacheIdentifier();
0524   unsigned long long cacheId3D = re3D.cacheIdentifier();
0525   if (cacheId2D != calibrationCacheId2D_ || cacheId3D != calibrationCacheId3D_)  //Calibration changed
0526   {
0527     ESHandle<TrackProbabilityCalibration> calib2DHandle = iSetup.getHandle(calib2d_token_);
0528     ESHandle<TrackProbabilityCalibration> calib3DHandle = iSetup.getHandle(calib3d_token_);
0529     probabilityEstimator_ =
0530         std::make_unique<HistogramProbabilityEstimator>(calib3DHandle.product(), calib2DHandle.product());
0531   }
0532 
0533   calibrationCacheId3D_ = cacheId3D;
0534   calibrationCacheId2D_ = cacheId2D;
0535 }
0536 
0537 //define this as a plug-in
0538 DEFINE_FWK_MODULE(DeepFlavourTagInfoProducer);