Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "CondFormats/DataRecord/interface/EcalPFSeedingThresholdsRcd.h"
0002 #include "CondFormats/EcalObjects/interface/EcalPFSeedingThresholds.h"
0003 #include "RecoParticleFlow/PFClusterProducer/interface/RecHitTopologicalCleanerBase.h"
0004 
0005 class ECALPFSeedCleaner : public RecHitTopologicalCleanerBase {
0006 public:
0007   ECALPFSeedCleaner(const edm::ParameterSet& conf, edm::ConsumesCollector& cc);
0008   ECALPFSeedCleaner(const ECALPFSeedCleaner&) = delete;
0009   ECALPFSeedCleaner& operator=(const ECALPFSeedCleaner&) = delete;
0010 
0011   void update(const edm::EventSetup&) override;
0012 
0013   void clean(const edm::Handle<reco::PFRecHitCollection>& input, std::vector<bool>& mask) override;
0014 
0015 private:
0016   edm::ESHandle<EcalPFSeedingThresholds> ths_;
0017   edm::ESGetToken<EcalPFSeedingThresholds, EcalPFSeedingThresholdsRcd> thsToken_;
0018 };
0019 
0020 DEFINE_EDM_PLUGIN(RecHitTopologicalCleanerFactory, ECALPFSeedCleaner, "ECALPFSeedCleaner");
0021 
0022 ECALPFSeedCleaner::ECALPFSeedCleaner(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0023     : RecHitTopologicalCleanerBase(conf, cc), thsToken_(cc.esConsumes<edm::Transition::BeginRun>()) {}
0024 
0025 void ECALPFSeedCleaner::update(const edm::EventSetup& iSetup) { ths_ = iSetup.getHandle(thsToken_); }
0026 
0027 void ECALPFSeedCleaner::clean(const edm::Handle<reco::PFRecHitCollection>& input, std::vector<bool>& mask) {
0028   //need to run over energy sorted rechits, as this is order used in seeding step
0029   // this can cause ambiguity, isn't it better to index by detid ?
0030   auto const& hits = *input;
0031   std::vector<unsigned> ordered_hits(hits.size());
0032   for (unsigned i = 0; i < hits.size(); ++i)
0033     ordered_hits[i] = i;
0034 
0035   std::sort(ordered_hits.begin(), ordered_hits.end(), [&](unsigned i, unsigned j) {
0036     return hits[i].energy() > hits[j].energy();
0037   });
0038 
0039   for (const auto& idx : ordered_hits) {
0040     if (!mask[idx])
0041       continue;  // is it useful ?
0042     const reco::PFRecHit& rechit = hits[idx];
0043 
0044     float threshold = (*ths_)[rechit.detId()];
0045     if (rechit.energy() < threshold)
0046       mask[idx] = false;
0047 
0048   }  // rechit loop
0049 }