Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-19 05:14:17

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/Utilities/interface/StreamID.h"
0009 #include "FWCore/Utilities/interface/ESGetToken.h"
0010 
0011 #include "DataFormats/PatCandidates/interface/Jet.h"
0012 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0013 
0014 #include "DataFormats/BTauReco/interface/ShallowTagInfo.h"
0015 
0016 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0017 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0018 
0019 #include "DataFormats/BTauReco/interface/UnifiedParticleTransformerAK4TagInfo.h"
0020 #include "DataFormats/BTauReco/interface/UnifiedParticleTransformerAK4Features.h"
0021 
0022 #include "RecoBTag/FeatureTools/interface/JetConverter.h"
0023 #include "RecoBTag/FeatureTools/interface/ShallowTagInfoConverter.h"
0024 #include "RecoBTag/FeatureTools/interface/SecondaryVertexConverter.h"
0025 #include "RecoBTag/FeatureTools/interface/NeutralCandidateConverter.h"
0026 #include "RecoBTag/FeatureTools/interface/ChargedCandidateConverter.h"
0027 #include "RecoBTag/FeatureTools/interface/LostTracksConverter.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 UnifiedParticleTransformerAK4TagInfoProducer : public edm::stream::EDProducer<> {
0063 public:
0064   explicit UnifiedParticleTransformerAK4TagInfoProducer(const edm::ParameterSet&);
0065   ~UnifiedParticleTransformerAK4TagInfoProducer() override = default;
0066 
0067   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0068 
0069 private:
0070   typedef std::vector<reco::UnifiedParticleTransformerAK4TagInfo> UnifiedParticleTransformerAK4TagInfoCollection;
0071   typedef reco::VertexCompositePtrCandidateCollection SVCollection;
0072   typedef reco::VertexCollection VertexCollection;
0073   typedef edm::AssociationMap<edm::OneToOne<reco::JetView, reco::JetView>> JetMatchMap;
0074 
0075   void beginStream(edm::StreamID) override {}
0076   void produce(edm::Event&, const edm::EventSetup&) override;
0077   void endStream() override {}
0078 
0079   const double jet_radius_;
0080   const double min_candidate_pt_;
0081   const bool flip_;
0082   const bool sort_cand_by_pt_;
0083   const bool fix_lt_sorting_;
0084 
0085   const edm::EDGetTokenT<edm::View<reco::Jet>> jet_token_;
0086   const edm::EDGetTokenT<VertexCollection> vtx_token_;
0087   const edm::EDGetTokenT<edm::View<reco::Candidate>> lt_token_;
0088   const edm::EDGetTokenT<SVCollection> sv_token_;
0089   edm::EDGetTokenT<JetMatchMap> unsubjet_map_token_;
0090   edm::EDGetTokenT<edm::ValueMap<float>> puppi_value_map_token_;
0091   edm::EDGetTokenT<edm::ValueMap<int>> pvasq_value_map_token_;
0092   edm::EDGetTokenT<edm::Association<VertexCollection>> pvas_token_;
0093   const edm::EDGetTokenT<edm::View<reco::Candidate>> candidateToken_;
0094   const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> track_builder_token_;
0095   bool use_puppi_value_map_;
0096   bool use_pvasq_value_map_;
0097   bool use_unsubjet_map_;
0098 
0099   const bool fallback_puppi_weight_;
0100   const bool fallback_vertex_association_;
0101 
0102   const bool is_weighted_jet_;
0103 
0104   const double min_jet_pt_;
0105   const double max_jet_eta_;
0106 };
0107 
0108 UnifiedParticleTransformerAK4TagInfoProducer::UnifiedParticleTransformerAK4TagInfoProducer(
0109     const edm::ParameterSet& iConfig)
0110     : jet_radius_(iConfig.getParameter<double>("jet_radius")),
0111       min_candidate_pt_(iConfig.getParameter<double>("min_candidate_pt")),
0112       flip_(iConfig.getParameter<bool>("flip")),
0113       sort_cand_by_pt_(iConfig.getParameter<bool>("sort_cand_by_pt")),
0114       fix_lt_sorting_(iConfig.getParameter<bool>("fix_lt_sorting")),
0115       jet_token_(consumes<edm::View<reco::Jet>>(iConfig.getParameter<edm::InputTag>("jets"))),
0116       vtx_token_(consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
0117       lt_token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("losttracks"))),
0118       sv_token_(consumes<SVCollection>(iConfig.getParameter<edm::InputTag>("secondary_vertices"))),
0119       candidateToken_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("candidates"))),
0120       track_builder_token_(
0121           esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"))),
0122       use_puppi_value_map_(false),
0123       use_pvasq_value_map_(false),
0124       use_unsubjet_map_(false),
0125       fallback_puppi_weight_(iConfig.getParameter<bool>("fallback_puppi_weight")),
0126       fallback_vertex_association_(iConfig.getParameter<bool>("fallback_vertex_association")),
0127       is_weighted_jet_(iConfig.getParameter<bool>("is_weighted_jet")),
0128       min_jet_pt_(iConfig.getParameter<double>("min_jet_pt")),
0129       max_jet_eta_(iConfig.getParameter<double>("max_jet_eta")) {
0130   produces<UnifiedParticleTransformerAK4TagInfoCollection>();
0131 
0132   const auto& puppi_value_map_tag = iConfig.getParameter<edm::InputTag>("puppi_value_map");
0133   if (!puppi_value_map_tag.label().empty()) {
0134     puppi_value_map_token_ = consumes<edm::ValueMap<float>>(puppi_value_map_tag);
0135     use_puppi_value_map_ = true;
0136   } else if (is_weighted_jet_) {
0137     throw edm::Exception(edm::errors::Configuration,
0138                          "puppi_value_map is not set but jet is weighted. Must set puppi_value_map.");
0139   }
0140 
0141   const auto& pvas_tag = iConfig.getParameter<edm::InputTag>("vertex_associator");
0142   if (!pvas_tag.label().empty()) {
0143     pvasq_value_map_token_ = consumes<edm::ValueMap<int>>(pvas_tag);
0144     pvas_token_ = consumes<edm::Association<VertexCollection>>(pvas_tag);
0145     use_pvasq_value_map_ = true;
0146   }
0147 
0148   const auto& unsubjet_map_tag = iConfig.getParameter<edm::InputTag>("unsubjet_map");
0149   if (!unsubjet_map_tag.label().empty()) {
0150     unsubjet_map_token_ = consumes<JetMatchMap>(unsubjet_map_tag);
0151     use_unsubjet_map_ = true;
0152   }
0153 }
0154 
0155 void UnifiedParticleTransformerAK4TagInfoProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0156   // pfUnifiedParticleTransformerAK4TagInfos
0157   edm::ParameterSetDescription desc;
0158   desc.add<double>("jet_radius", 0.4);
0159   desc.add<double>("min_candidate_pt", 0.10);
0160   desc.add<bool>("flip", false);
0161   desc.add<bool>("sort_cand_by_pt", false);
0162   desc.add<bool>("fix_lt_sorting", false);
0163   desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
0164   desc.add<edm::InputTag>("losttracks", edm::InputTag("lostTracks"));
0165   desc.add<edm::InputTag>("puppi_value_map", edm::InputTag("puppi"));
0166   desc.add<edm::InputTag>("secondary_vertices", edm::InputTag("inclusiveCandidateSecondaryVertices"));
0167   desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
0168   desc.add<edm::InputTag>("unsubjet_map", {});
0169   desc.add<edm::InputTag>("candidates", edm::InputTag("packedPFCandidates"));
0170   desc.add<edm::InputTag>("vertex_associator", edm::InputTag("primaryVertexAssociation", "original"));
0171   desc.add<bool>("fallback_puppi_weight", false);
0172   desc.add<bool>("fallback_vertex_association", false);
0173   desc.add<bool>("is_weighted_jet", false);
0174   desc.add<double>("min_jet_pt", 0.0);
0175   desc.add<double>("max_jet_eta", 2.5);
0176   descriptions.add("pfUnifiedParticleTransformerAK4TagInfos", desc);
0177 }
0178 
0179 void UnifiedParticleTransformerAK4TagInfoProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0180   auto output_tag_infos = std::make_unique<UnifiedParticleTransformerAK4TagInfoCollection>();
0181   edm::Handle<edm::View<reco::Jet>> jets;
0182   iEvent.getByToken(jet_token_, jets);
0183 
0184   edm::Handle<JetMatchMap> unsubjet_map;
0185   if (use_unsubjet_map_)
0186     iEvent.getByToken(unsubjet_map_token_, unsubjet_map);
0187 
0188   edm::Handle<VertexCollection> vtxs;
0189   iEvent.getByToken(vtx_token_, vtxs);
0190   if (vtxs->empty()) {
0191     // produce empty TagInfos in case no primary vertex
0192     iEvent.put(std::move(output_tag_infos));
0193     return;  // exit event
0194   }
0195 
0196   edm::Handle<edm::View<reco::Candidate>> LTs;
0197   iEvent.getByToken(lt_token_, LTs);
0198 
0199   // reference to primary vertex
0200   const auto& pv = vtxs->at(0);
0201 
0202   edm::Handle<edm::View<reco::Candidate>> tracks;
0203   iEvent.getByToken(candidateToken_, tracks);
0204 
0205   edm::Handle<SVCollection> svs;
0206   iEvent.getByToken(sv_token_, svs);
0207 
0208   edm::Handle<edm::ValueMap<float>> puppi_value_map;
0209   if (use_puppi_value_map_) {
0210     iEvent.getByToken(puppi_value_map_token_, puppi_value_map);
0211   }
0212 
0213   edm::Handle<edm::ValueMap<int>> pvasq_value_map;
0214   edm::Handle<edm::Association<VertexCollection>> pvas;
0215   if (use_pvasq_value_map_) {
0216     iEvent.getByToken(pvasq_value_map_token_, pvasq_value_map);
0217     iEvent.getByToken(pvas_token_, pvas);
0218   }
0219 
0220   edm::ESHandle<TransientTrackBuilder> track_builder = iSetup.getHandle(track_builder_token_);
0221 
0222   for (std::size_t jet_n = 0; jet_n < jets->size(); jet_n++) {
0223     // create data containing structure
0224     btagbtvdeep::UnifiedParticleTransformerAK4Features features;
0225 
0226     // reco jet reference (use as much as possible)
0227     const auto& jet = jets->at(jet_n);
0228     if (jet.pt() < min_jet_pt_) {
0229       features.is_filled = false;
0230     }
0231     if (std::abs(jet.eta()) > max_jet_eta_) {
0232       features.is_filled = false;
0233     }
0234     // dynamical casting to pointers, null if not possible
0235     const auto* pf_jet = dynamic_cast<const reco::PFJet*>(&jet);
0236     const auto* pat_jet = dynamic_cast<const pat::Jet*>(&jet);
0237     edm::RefToBase<reco::Jet> jet_ref(jets, jet_n);
0238     const auto& unsubJet =
0239         (use_unsubjet_map_ && (*unsubjet_map)[jet_ref].isNonnull()) ? *(*unsubjet_map)[jet_ref] : jet;
0240 
0241     if (features.is_filled) {
0242       math::XYZVector jet_dir = jet.momentum().Unit();
0243       GlobalVector jet_ref_track_dir(jet.px(), jet.py(), jet.pz());
0244 
0245       // copy which will be sorted
0246       auto svs_sorted = *svs;
0247       // sort by dxy
0248       std::sort(svs_sorted.begin(), svs_sorted.end(), [&pv](const auto& sva, const auto& svb) {
0249         return btagbtvdeep::sv_vertex_comparator(sva, svb, pv);
0250       });
0251       // fill features from secondary vertices
0252       for (const auto& sv : svs_sorted) {
0253         if (reco::deltaR2(sv, jet_dir) > (jet_radius_ * jet_radius_))
0254           continue;
0255         else {
0256           features.sv_features.emplace_back();
0257           // in C++17 could just get from emplace_back output
0258           auto& sv_features = features.sv_features.back();
0259           btagbtvdeep::svToFeatures(sv, pv, jet, sv_features, flip_);
0260         }
0261       }
0262 
0263       // stuff required for dealing with pf candidates and lost tracks
0264       std::vector<btagbtvdeep::SortingClass<size_t>> c_sorted, n_sorted, lt_sorted;
0265       // to cache the TrackInfo
0266       std::map<unsigned int, btagbtvdeep::TrackInfoBuilder> trackinfos;
0267       std::map<unsigned int, btagbtvdeep::TrackInfoBuilder> lt_trackinfos;
0268       // unsorted reference to sv
0269       const auto& svs_unsorted = *svs;
0270       std::vector<reco::CandidatePtr> ltPtrs;
0271       static constexpr size_t max_lt_ = 5;
0272 
0273       //Adding the lost tracks associated with the jets
0274       for (size_t i = 0; i < LTs->size(); ++i) {
0275         auto cand = LTs->ptrAt(i);
0276         if ((reco::deltaR(*cand, jet) < 0.2)) {
0277           const auto* PackedCandidate_ = dynamic_cast<const pat::PackedCandidate*>(&(*cand));
0278           if (PackedCandidate_) {
0279             if (PackedCandidate_->pt() < 1.0)
0280               continue;
0281             auto& trackinfo = lt_trackinfos.emplace(i, track_builder).first->second;
0282             trackinfo.buildTrackInfo(PackedCandidate_, jet_dir, jet_ref_track_dir, pv);
0283 
0284             if (sort_cand_by_pt_)
0285               lt_sorted.emplace_back(i,
0286                                      PackedCandidate_->pt() / jet.pt(),
0287                                      trackinfo.getTrackSip2dSig(),
0288                                      -btagbtvdeep::mindrsvpfcand(svs_unsorted, PackedCandidate_));
0289             else
0290               lt_sorted.emplace_back(i,
0291                                      trackinfo.getTrackSip2dSig(),
0292                                      -btagbtvdeep::mindrsvpfcand(svs_unsorted, PackedCandidate_),
0293                                      PackedCandidate_->pt() / jet.pt());
0294 
0295             ltPtrs.push_back(cand);
0296           }
0297         }
0298       }
0299 
0300       // sort lt collection
0301       std::sort(lt_sorted.begin(), lt_sorted.end(), btagbtvdeep::SortingClass<std::size_t>::compareByABCInv);
0302       int n_ltcand_ = std::min(lt_sorted.size(), max_lt_);
0303 
0304       std::vector<size_t> lt_sortedindices;
0305       lt_sortedindices = btagbtvdeep::invertSortingVector(lt_sorted);
0306 
0307       // set right size to vectors
0308       features.lt_features.clear();
0309       features.lt_features.resize(lt_sorted.size());
0310 
0311       for (unsigned int i = 0; i < (unsigned int)n_ltcand_; i++) {
0312         //auto cand = LTs->ptrAt(i);
0313         //const auto *PackedCandidate_ = dynamic_cast<const pat::PackedCandidate*>(&(*cand));
0314         const auto* PackedCandidate_ = dynamic_cast<const pat::PackedCandidate*>(&(*ltPtrs.at(i)));
0315         if (!PackedCandidate_)
0316           continue;
0317         if (PackedCandidate_->pt() < min_candidate_pt_)
0318           continue;
0319 
0320         if (PackedCandidate_->charge() != 0) {
0321           //auto reco_cand = dynamic_cast<const reco::PFCandidate*>(cand);
0322           float puppiw = PackedCandidate_->puppiWeight();
0323 
0324           float drminpfcandsv = btagbtvdeep::mindrsvpfcand(svs_unsorted, PackedCandidate_);
0325           float distminpfcandsv = 0;
0326 
0327           size_t entry = lt_sortedindices.at(fix_lt_sorting_ ? lt_sorted[i].get() : i);
0328           // get cached track info
0329           auto& trackinfo = lt_trackinfos.emplace(i, track_builder).first->second;
0330           trackinfo.buildTrackInfo(PackedCandidate_, jet_dir, jet_ref_track_dir, pv);
0331           // get_ref to vector element
0332           auto& lt_features = features.lt_features.at(entry);
0333 
0334           if (PackedCandidate_) {
0335             if (PackedCandidate_->hasTrackDetails()) {
0336               const reco::Track& PseudoTrack = PackedCandidate_->pseudoTrack();
0337               reco::TransientTrack transientTrack;
0338               transientTrack = track_builder->build(PseudoTrack);
0339               distminpfcandsv = btagbtvdeep::mindistsvpfcand(svs_unsorted, transientTrack);
0340             }
0341 
0342             btagbtvdeep::packedCandidateToFeatures(PackedCandidate_,
0343                                                    jet,
0344                                                    trackinfo,
0345                                                    is_weighted_jet_,
0346                                                    drminpfcandsv,
0347                                                    static_cast<float>(jet_radius_),
0348                                                    puppiw,
0349                                                    lt_features,
0350                                                    flip_,
0351                                                    distminpfcandsv);
0352           }
0353         }
0354       }
0355 
0356       // fill collection, from DeepTNtuples plus some styling
0357       for (unsigned int i = 0; i < unsubJet.numberOfDaughters(); i++) {
0358         auto cand = unsubJet.daughter(i);
0359         if (cand) {
0360           // candidates under 100MeV (configurable) are not considered
0361           // might change if we use also white-listing
0362           if (cand->pt() < min_candidate_pt_)
0363             continue;
0364           if (cand->charge() != 0) {
0365             auto& trackinfo = trackinfos.emplace(i, track_builder).first->second;
0366             trackinfo.buildTrackInfo(cand, jet_dir, jet_ref_track_dir, pv);
0367 
0368             if (sort_cand_by_pt_)
0369               c_sorted.emplace_back(i,
0370                                     cand->pt() / jet.pt(),
0371                                     trackinfo.getTrackSip2dSig(),
0372                                     -btagbtvdeep::mindrsvpfcand(svs_unsorted, cand));
0373             else
0374               c_sorted.emplace_back(i,
0375                                     trackinfo.getTrackSip2dSig(),
0376                                     -btagbtvdeep::mindrsvpfcand(svs_unsorted, cand),
0377                                     cand->pt() / jet.pt());
0378           } else {
0379             if (sort_cand_by_pt_)
0380               n_sorted.emplace_back(i, cand->pt() / jet.pt(), -btagbtvdeep::mindrsvpfcand(svs_unsorted, cand), -1);
0381             else
0382               n_sorted.emplace_back(i, -1, -btagbtvdeep::mindrsvpfcand(svs_unsorted, cand), cand->pt() / jet.pt());
0383           }
0384         }
0385       }
0386 
0387       // sort collections (open the black-box if you please)
0388       std::sort(c_sorted.begin(), c_sorted.end(), btagbtvdeep::SortingClass<std::size_t>::compareByABCInv);
0389       std::sort(n_sorted.begin(), n_sorted.end(), btagbtvdeep::SortingClass<std::size_t>::compareByABCInv);
0390 
0391       std::vector<size_t> c_sortedindices, n_sortedindices;
0392 
0393       // this puts 0 everywhere and the right position in ind
0394       c_sortedindices = btagbtvdeep::invertSortingVector(c_sorted);
0395       n_sortedindices = btagbtvdeep::invertSortingVector(n_sorted);
0396 
0397       // set right size to vectors
0398       features.c_pf_features.clear();
0399       features.c_pf_features.resize(c_sorted.size());
0400       features.n_pf_features.clear();
0401       features.n_pf_features.resize(n_sorted.size());
0402 
0403       for (unsigned int i = 0; i < unsubJet.numberOfDaughters(); i++) {
0404         // get pointer and check that is correct
0405         auto cand = dynamic_cast<const reco::Candidate*>(unsubJet.daughter(i));
0406         if (!cand)
0407           continue;
0408         // candidates under 100MeV are not considered
0409         // might change if we use also white-listing
0410         if (cand->pt() < 0.100)
0411           continue;
0412 
0413         auto packed_cand = dynamic_cast<const pat::PackedCandidate*>(cand);
0414         auto reco_cand = dynamic_cast<const reco::PFCandidate*>(cand);
0415 
0416         // need some edm::Ptr or edm::Ref if reco candidates
0417         reco::PFCandidatePtr reco_ptr;
0418         if (pf_jet) {
0419           reco_ptr = pf_jet->getPFConstituent(i);
0420         } else if (pat_jet && reco_cand) {
0421           reco_ptr = pat_jet->getPFConstituent(i);
0422         }
0423 
0424         reco::CandidatePtr cand_ptr;
0425         if (pat_jet) {
0426           cand_ptr = pat_jet->sourceCandidatePtr(i);
0427         }
0428 
0429         //
0430         // Access puppi weight from ValueMap.
0431         //
0432         float puppiw = 1.0;  // Set to fallback value
0433 
0434         if (reco_cand) {
0435           if (use_puppi_value_map_)
0436             puppiw = (*puppi_value_map)[reco_ptr];
0437           else if (!fallback_puppi_weight_) {
0438             throw edm::Exception(edm::errors::InvalidReference, "PUPPI value map missing")
0439                 << "use fallback_puppi_weight option to use " << puppiw << " for reco_cand as default";
0440           }
0441         } else if (packed_cand) {
0442           if (use_puppi_value_map_)
0443             puppiw = (*puppi_value_map)[cand_ptr];
0444           else if (!fallback_puppi_weight_) {
0445             throw edm::Exception(edm::errors::InvalidReference, "PUPPI value map missing")
0446                 << "use fallback_puppi_weight option to use " << puppiw << " for packed_cand as default";
0447           }
0448         } else {
0449           throw edm::Exception(edm::errors::InvalidReference)
0450               << "Cannot convert to either reco::PFCandidate or pat::PackedCandidate";
0451         }
0452 
0453         float drminpfcandsv = btagbtvdeep::mindrsvpfcand(svs_unsorted, cand);
0454         float distminpfcandsv = 0;
0455         if (cand->charge() != 0) {
0456           // is charged candidate
0457           auto entry = c_sortedindices.at(i);
0458 
0459           // get cached track info
0460           auto& trackinfo = trackinfos.at(i);
0461 
0462           // get_ref to vector element
0463           auto& c_pf_features = features.c_pf_features.at(entry);
0464           // fill feature structure
0465           if (packed_cand) {
0466             if (packed_cand->hasTrackDetails()) {
0467               const reco::Track& PseudoTrack = packed_cand->pseudoTrack();
0468               reco::TransientTrack transientTrack;
0469               transientTrack = track_builder->build(PseudoTrack);
0470               distminpfcandsv = btagbtvdeep::mindistsvpfcand(svs_unsorted, transientTrack);
0471             }
0472 
0473             btagbtvdeep::packedCandidateToFeatures(packed_cand,
0474                                                    jet,
0475                                                    trackinfo,
0476                                                    is_weighted_jet_,
0477                                                    drminpfcandsv,
0478                                                    static_cast<float>(jet_radius_),
0479                                                    puppiw,
0480                                                    c_pf_features,
0481                                                    flip_,
0482                                                    distminpfcandsv);
0483           } else if (reco_cand) {
0484             // get vertex association quality
0485             int pv_ass_quality = 0;  // fallback value
0486             if (use_pvasq_value_map_) {
0487               pv_ass_quality = (*pvasq_value_map)[reco_ptr];
0488             } else if (!fallback_vertex_association_) {
0489               throw edm::Exception(edm::errors::InvalidReference, "vertex association missing")
0490                   << "use fallback_vertex_association option to use" << pv_ass_quality
0491                   << "as default quality and closest dz PV as criteria";
0492             }
0493             // getting the PV as PackedCandidatesProducer
0494             // but using not the slimmed but original vertices
0495             auto ctrack = reco_cand->bestTrack();
0496             int pvi = -1;
0497             float dist = 1e99;
0498             for (size_t ii = 0; ii < vtxs->size(); ii++) {
0499               float dz = (ctrack) ? std::abs(ctrack->dz(((*vtxs)[ii]).position())) : 0;
0500               if (dz < dist) {
0501                 pvi = ii;
0502                 dist = dz;
0503               }
0504             }
0505             auto PV = reco::VertexRef(vtxs, pvi);
0506             if (use_pvasq_value_map_) {
0507               const reco::VertexRef& PV_orig = (*pvas)[reco_ptr];
0508               if (PV_orig.isNonnull())
0509                 PV = reco::VertexRef(vtxs, PV_orig.key());
0510             }
0511             btagbtvdeep::recoCandidateToFeatures(reco_cand,
0512                                                  jet,
0513                                                  trackinfo,
0514                                                  is_weighted_jet_,
0515                                                  drminpfcandsv,
0516                                                  static_cast<float>(jet_radius_),
0517                                                  puppiw,
0518                                                  pv_ass_quality,
0519                                                  PV,
0520                                                  c_pf_features,
0521                                                  flip_,
0522                                                  distminpfcandsv);
0523           }
0524         } else {
0525           // is neutral candidate
0526           auto entry = n_sortedindices.at(i);
0527           // get_ref to vector element
0528           auto& n_pf_features = features.n_pf_features.at(entry);
0529           // fill feature structure
0530           if (packed_cand) {
0531             btagbtvdeep::packedCandidateToFeatures(packed_cand,
0532                                                    jet,
0533                                                    is_weighted_jet_,
0534                                                    drminpfcandsv,
0535                                                    static_cast<float>(jet_radius_),
0536                                                    puppiw,
0537                                                    n_pf_features);
0538           } else if (reco_cand) {
0539             btagbtvdeep::recoCandidateToFeatures(
0540                 reco_cand, jet, is_weighted_jet_, drminpfcandsv, static_cast<float>(jet_radius_), puppiw, n_pf_features);
0541           }
0542         }
0543       }
0544     }
0545     output_tag_infos->emplace_back(features, jet_ref);
0546   }
0547   iEvent.put(std::move(output_tag_infos));
0548 }
0549 
0550 //define this as a plug-in
0551 DEFINE_FWK_MODULE(UnifiedParticleTransformerAK4TagInfoProducer);