Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // system include files
0002 #include <memory>
0003 
0004 // user include files
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/stream/EDProducer.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0013 #include "DataFormats/Scouting/interface/Run3ScoutingParticle.h"
0014 #include "DataFormats/Common/interface/OrphanHandle.h"
0015 
0016 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
0017 #include "fastjet/contrib/SoftKiller.hh"
0018 
0019 class Run3ScoutingParticleToRecoPFCandidateProducer : public edm::stream::EDProducer<> {
0020 public:
0021   explicit Run3ScoutingParticleToRecoPFCandidateProducer(const edm::ParameterSet &);
0022   ~Run3ScoutingParticleToRecoPFCandidateProducer() override;
0023 
0024   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0025   void beginStream(edm::StreamID) override {}
0026   void produce(edm::Event &iEvent, edm::EventSetup const &setup) override;
0027   void endStream() override {}
0028 
0029   void createPFCandidates(edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
0030                           std::unique_ptr<reco::PFCandidateCollection> &pfcands);
0031   void createPFCandidatesSK(edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
0032                             std::unique_ptr<reco::PFCandidateCollection> &pfcands);
0033   reco::PFCandidate createPFCand(Run3ScoutingParticle scoutingparticle);
0034   void clearVars();
0035 
0036 private:
0037   const edm::EDGetTokenT<std::vector<Run3ScoutingParticle>> input_scoutingparticle_token_;
0038   const edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> particletable_token_;
0039   bool use_softKiller_;
0040   bool use_CHS_;
0041   const HepPDT::ParticleDataTable *pdTable_;
0042 
0043   std::vector<int8_t> vertexIndex_;
0044   std::vector<float> normchi2_;
0045   std::vector<float> dz_;
0046   std::vector<float> dxy_;
0047   std::vector<float> dzsig_;
0048   std::vector<float> dxysig_;
0049   std::vector<int> lostInnerHits_;
0050   std::vector<int> quality_;
0051   std::vector<float> trkPt_;
0052   std::vector<float> trkEta_;
0053   std::vector<float> trkPhi_;
0054 };
0055 
0056 //
0057 // constructors and destructor
0058 //
0059 Run3ScoutingParticleToRecoPFCandidateProducer::Run3ScoutingParticleToRecoPFCandidateProducer(
0060     const edm::ParameterSet &iConfig)
0061     : input_scoutingparticle_token_(consumes(iConfig.getParameter<edm::InputTag>("scoutingparticle"))),
0062       particletable_token_(esConsumes<HepPDT::ParticleDataTable, edm::DefaultRecord>()),
0063       use_softKiller_(iConfig.getParameter<bool>("softKiller")),
0064       use_CHS_(iConfig.getParameter<bool>("CHS")) {
0065   //register products
0066   produces<reco::PFCandidateCollection>();
0067   produces<edm::ValueMap<int>>("vertexIndex");
0068   produces<edm::ValueMap<float>>("normchi2");
0069   produces<edm::ValueMap<float>>("dz");
0070   produces<edm::ValueMap<float>>("dxy");
0071   produces<edm::ValueMap<float>>("dzsig");
0072   produces<edm::ValueMap<float>>("dxysig");
0073   produces<edm::ValueMap<int>>("lostInnerHits");
0074   produces<edm::ValueMap<int>>("quality");
0075   produces<edm::ValueMap<float>>("trkPt");
0076   produces<edm::ValueMap<float>>("trkEta");
0077   produces<edm::ValueMap<float>>("trkPhi");
0078 }
0079 
0080 Run3ScoutingParticleToRecoPFCandidateProducer::~Run3ScoutingParticleToRecoPFCandidateProducer() = default;
0081 
0082 reco::PFCandidate Run3ScoutingParticleToRecoPFCandidateProducer::createPFCand(Run3ScoutingParticle scoutingparticle) {
0083   auto m = pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId())) != nullptr
0084                ? pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId()))->mass()
0085                : -99.f;
0086   auto q = pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId())) != nullptr
0087                ? pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId()))->charge()
0088                : -99.f;
0089   if (m < -90 or q < -90) {
0090     LogDebug("createPFCand") << "<Run3ScoutingParticleToRecoPFCandidateProducer::createPFCand>:" << std::endl
0091                              << "Unrecognisable pdgId - skipping particle" << std::endl;
0092     return reco::PFCandidate();
0093   }
0094 
0095   float px = scoutingparticle.pt() * cos(scoutingparticle.phi());
0096   float py = scoutingparticle.pt() * sin(scoutingparticle.phi());
0097   float pz = scoutingparticle.pt() * sinh(scoutingparticle.eta());
0098   float p = scoutingparticle.pt() * cosh(scoutingparticle.eta());
0099   float energy = std::sqrt(p * p + m * m);
0100   reco::Particle::LorentzVector p4(px, py, pz, energy);
0101 
0102   static const reco::PFCandidate dummy;
0103   auto pfcand = reco::PFCandidate(q, p4, dummy.translatePdgIdToType(scoutingparticle.pdgId()));
0104 
0105   bool relativeTrackVars = scoutingparticle.relative_trk_vars();
0106   vertexIndex_.push_back(scoutingparticle.vertex());
0107   normchi2_.push_back(scoutingparticle.normchi2());
0108   dz_.push_back(scoutingparticle.dz());
0109   dxy_.push_back(scoutingparticle.dxy());
0110   dzsig_.push_back(scoutingparticle.dzsig());
0111   dxysig_.push_back(scoutingparticle.dxysig());
0112   lostInnerHits_.push_back(scoutingparticle.lostInnerHits());
0113   quality_.push_back(scoutingparticle.quality());
0114   trkPt_.push_back(relativeTrackVars ? scoutingparticle.trk_pt() + scoutingparticle.pt() : scoutingparticle.trk_pt());
0115   trkEta_.push_back(relativeTrackVars ? scoutingparticle.trk_eta() + scoutingparticle.eta()
0116                                       : scoutingparticle.trk_eta());
0117   trkPhi_.push_back(relativeTrackVars ? scoutingparticle.trk_phi() + scoutingparticle.phi()
0118                                       : scoutingparticle.trk_phi());
0119 
0120   return pfcand;
0121 }
0122 
0123 void Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidates(
0124     edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
0125     std::unique_ptr<reco::PFCandidateCollection> &pfcands) {
0126   for (unsigned int icand = 0; icand < scoutingparticleHandle->size(); ++icand) {
0127     auto &scoutingparticle = (*scoutingparticleHandle)[icand];
0128 
0129     if (use_CHS_ and scoutingparticle.vertex() > 0)
0130       continue;
0131 
0132     auto pfcand = createPFCand(scoutingparticle);
0133     if (pfcand.energy() != 0)
0134       pfcands->push_back(pfcand);
0135   }
0136 }
0137 
0138 void Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidatesSK(
0139     edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
0140     std::unique_ptr<reco::PFCandidateCollection> &pfcands) {
0141   std::vector<fastjet::PseudoJet> fj;
0142 
0143   for (auto iter = scoutingparticleHandle->begin(),
0144             ibegin = scoutingparticleHandle->begin(),
0145             iend = scoutingparticleHandle->end();
0146        iter != iend;
0147        ++iter) {
0148     auto m = pdTable_->particle(HepPDT::ParticleID(iter->pdgId())) != nullptr
0149                  ? pdTable_->particle(HepPDT::ParticleID(iter->pdgId()))->mass()
0150                  : -99.f;
0151     if (m < -90) {
0152       LogDebug("createPFCandidatesSK") << "<Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidatesSK>:"
0153                                        << std::endl
0154                                        << "Unrecognisable pdgId - skipping particle" << std::endl;
0155       continue;
0156     }
0157     math::PtEtaPhiMLorentzVector p4(iter->pt(), iter->eta(), iter->phi(), m);
0158     fj.push_back(fastjet::PseudoJet(p4.px(), p4.py(), p4.pz(), p4.energy()));
0159     fj.back().set_user_index(iter - ibegin);
0160   }
0161 
0162   fastjet::contrib::SoftKiller soft_killer(5, 0.4);
0163   std::vector<fastjet::PseudoJet> soft_killed_particles = soft_killer(fj);
0164 
0165   for (auto &particle : soft_killed_particles) {
0166     const Run3ScoutingParticle scoutingparticle = scoutingparticleHandle->at(particle.user_index());
0167     auto pfcand = createPFCand(scoutingparticle);
0168     if (pfcand.energy() != 0)
0169       pfcands->push_back(pfcand);
0170   }
0171 }
0172 
0173 // ------------ method called to produce the data  ------------
0174 void Run3ScoutingParticleToRecoPFCandidateProducer::produce(edm::Event &iEvent, edm::EventSetup const &setup) {
0175   using namespace edm;
0176 
0177   auto pdt = setup.getHandle(particletable_token_);
0178   pdTable_ = pdt.product();
0179 
0180   Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle;
0181   iEvent.getByToken(input_scoutingparticle_token_, scoutingparticleHandle);
0182 
0183   auto pfcands = std::make_unique<reco::PFCandidateCollection>();
0184 
0185   if (use_softKiller_) {
0186     createPFCandidatesSK(scoutingparticleHandle, pfcands);
0187   } else {
0188     createPFCandidates(scoutingparticleHandle, pfcands);
0189   }
0190 
0191   edm::OrphanHandle<reco::PFCandidateCollection> oh = iEvent.put(std::move(pfcands));
0192 
0193   std::unique_ptr<edm::ValueMap<int>> vertexIndex_VM(new edm::ValueMap<int>());
0194   edm::ValueMap<int>::Filler filler_vertexIndex(*vertexIndex_VM);
0195   filler_vertexIndex.insert(oh, vertexIndex_.begin(), vertexIndex_.end());
0196   filler_vertexIndex.fill();
0197   iEvent.put(std::move(vertexIndex_VM), "vertexIndex");
0198 
0199   std::unique_ptr<edm::ValueMap<float>> normchi2_VM(new edm::ValueMap<float>());
0200   edm::ValueMap<float>::Filler filler_normchi2(*normchi2_VM);
0201   filler_normchi2.insert(oh, normchi2_.begin(), normchi2_.end());
0202   filler_normchi2.fill();
0203   iEvent.put(std::move(normchi2_VM), "normchi2");
0204 
0205   std::unique_ptr<edm::ValueMap<float>> dz_VM(new edm::ValueMap<float>());
0206   edm::ValueMap<float>::Filler filler_dz(*dz_VM);
0207   filler_dz.insert(oh, dz_.begin(), dz_.end());
0208   filler_dz.fill();
0209   iEvent.put(std::move(dz_VM), "dz");
0210 
0211   std::unique_ptr<edm::ValueMap<float>> dxy_VM(new edm::ValueMap<float>());
0212   edm::ValueMap<float>::Filler filler_dxy(*dxy_VM);
0213   filler_dxy.insert(oh, dxy_.begin(), dxy_.end());
0214   filler_dxy.fill();
0215   iEvent.put(std::move(dxy_VM), "dxy");
0216 
0217   std::unique_ptr<edm::ValueMap<float>> dzsig_VM(new edm::ValueMap<float>());
0218   edm::ValueMap<float>::Filler filler_dzsig(*dzsig_VM);
0219   filler_dzsig.insert(oh, dzsig_.begin(), dzsig_.end());
0220   filler_dzsig.fill();
0221   iEvent.put(std::move(dzsig_VM), "dzsig");
0222 
0223   std::unique_ptr<edm::ValueMap<float>> dxysig_VM(new edm::ValueMap<float>());
0224   edm::ValueMap<float>::Filler filler_dxysig(*dxysig_VM);
0225   filler_dxysig.insert(oh, dxysig_.begin(), dxysig_.end());
0226   filler_dxysig.fill();
0227   iEvent.put(std::move(dxysig_VM), "dxysig");
0228 
0229   std::unique_ptr<edm::ValueMap<int>> lostInnerHits_VM(new edm::ValueMap<int>());
0230   edm::ValueMap<int>::Filler filler_lostInnerHits(*lostInnerHits_VM);
0231   filler_lostInnerHits.insert(oh, lostInnerHits_.begin(), lostInnerHits_.end());
0232   filler_lostInnerHits.fill();
0233   iEvent.put(std::move(lostInnerHits_VM), "lostInnerHits");
0234 
0235   std::unique_ptr<edm::ValueMap<int>> quality_VM(new edm::ValueMap<int>());
0236   edm::ValueMap<int>::Filler filler_quality(*quality_VM);
0237   filler_quality.insert(oh, quality_.begin(), quality_.end());
0238   filler_quality.fill();
0239   iEvent.put(std::move(quality_VM), "quality");
0240 
0241   std::unique_ptr<edm::ValueMap<float>> trkPt_VM(new edm::ValueMap<float>());
0242   edm::ValueMap<float>::Filler filler_trkPt(*trkPt_VM);
0243   filler_trkPt.insert(oh, trkPt_.begin(), trkPt_.end());
0244   filler_trkPt.fill();
0245   iEvent.put(std::move(trkPt_VM), "trkPt");
0246 
0247   std::unique_ptr<edm::ValueMap<float>> trkEta_VM(new edm::ValueMap<float>());
0248   edm::ValueMap<float>::Filler filler_trkEta(*trkEta_VM);
0249   filler_trkEta.insert(oh, trkEta_.begin(), trkEta_.end());
0250   filler_trkEta.fill();
0251   iEvent.put(std::move(trkEta_VM), "trkEta");
0252 
0253   std::unique_ptr<edm::ValueMap<float>> trkPhi_VM(new edm::ValueMap<float>());
0254   edm::ValueMap<float>::Filler filler_trkPhi(*trkPhi_VM);
0255   filler_trkPhi.insert(oh, trkPhi_.begin(), trkPhi_.end());
0256   filler_trkPhi.fill();
0257   iEvent.put(std::move(trkPhi_VM), "trkPhi");
0258 
0259   clearVars();
0260 }
0261 
0262 void Run3ScoutingParticleToRecoPFCandidateProducer::clearVars() {
0263   vertexIndex_.clear();
0264   normchi2_.clear();
0265   dz_.clear();
0266   dxy_.clear();
0267   dzsig_.clear();
0268   dxysig_.clear();
0269   lostInnerHits_.clear();
0270   quality_.clear();
0271   trkPt_.clear();
0272   trkEta_.clear();
0273   trkPhi_.clear();
0274 }
0275 
0276 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0277 void Run3ScoutingParticleToRecoPFCandidateProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0278   edm::ParameterSetDescription desc;
0279   desc.add<edm::InputTag>("scoutingparticle", edm::InputTag("hltScoutingPFPacker"));
0280   desc.add<bool>("softKiller", false);
0281   desc.add<bool>("CHS", false);
0282   descriptions.addWithDefaultLabel(desc);
0283 }
0284 
0285 // declare this class as a framework plugin
0286 DEFINE_FWK_MODULE(Run3ScoutingParticleToRecoPFCandidateProducer);