Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DataFormats/ForwardDetId/interface/HGCalDetId.h"
0002 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFraction.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "RecoParticleFlow/PFClusterProducer/interface/InitialClusteringStepBase.h"
0006 #include "SimDataFormats/CaloAnalysis/interface/SimCluster.h"
0007 #include "SimDataFormats/CaloAnalysis/interface/SimClusterFwd.h"
0008 #include "CondFormats/DataRecord/interface/HcalPFCutsRcd.h"
0009 #include "CondTools/Hcal/interface/HcalPFCutsHandler.h"
0010 
0011 class GenericSimClusterMapper : public InitialClusteringStepBase {
0012   typedef GenericSimClusterMapper B2DGT;
0013 
0014 public:
0015   GenericSimClusterMapper(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0016       : InitialClusteringStepBase(conf, cc) {
0017     _simClusterToken = cc.consumes<SimClusterCollection>(conf.getParameter<edm::InputTag>("simClusterSrc"));
0018   }
0019   ~GenericSimClusterMapper() override = default;
0020   GenericSimClusterMapper(const B2DGT&) = delete;
0021   B2DGT& operator=(const B2DGT&) = delete;
0022 
0023   void updateEvent(const edm::Event&) final;
0024 
0025   void buildClusters(const edm::Handle<reco::PFRecHitCollection>&,
0026                      const std::vector<bool>&,
0027                      const std::vector<bool>&,
0028                      reco::PFClusterCollection&,
0029                      const HcalPFCuts*) override;
0030 
0031 private:
0032   edm::EDGetTokenT<SimClusterCollection> _simClusterToken;
0033   edm::Handle<SimClusterCollection> _simClusterH;
0034 };
0035 
0036 DEFINE_EDM_PLUGIN(InitialClusteringStepFactory, GenericSimClusterMapper, "GenericSimClusterMapper");
0037 
0038 #ifdef PFLOW_DEBUG
0039 #define LOGVERB(x) edm::LogVerbatim(x)
0040 #define LOGWARN(x) edm::LogWarning(x)
0041 #define LOGERR(x) edm::LogError(x)
0042 #define LOGDRESSED(x) edm::LogInfo(x)
0043 #else
0044 #define LOGVERB(x) LogTrace(x)
0045 #define LOGWARN(x) edm::LogWarning(x)
0046 #define LOGERR(x) edm::LogError(x)
0047 #define LOGDRESSED(x) LogDebug(x)
0048 #endif
0049 
0050 void GenericSimClusterMapper::updateEvent(const edm::Event& ev) { ev.getByToken(_simClusterToken, _simClusterH); }
0051 
0052 void GenericSimClusterMapper::buildClusters(const edm::Handle<reco::PFRecHitCollection>& input,
0053                                             const std::vector<bool>& rechitMask,
0054                                             const std::vector<bool>& seedable,
0055                                             reco::PFClusterCollection& output,
0056                                             const HcalPFCuts* hcalCuts) {
0057   const SimClusterCollection& simClusters = *_simClusterH;
0058   auto const& hits = *input;
0059 
0060   // for quick indexing back to hit energy
0061   std::unordered_map<uint32_t, size_t> detIdToIndex(hits.size());
0062   for (uint32_t i = 0; i < hits.size(); ++i) {
0063     detIdToIndex[hits[i].detId()] = i;
0064     auto ref = makeRefhit(input, i);
0065   }
0066 
0067   for (const auto& sc : simClusters) {
0068     output.emplace_back();
0069     reco::PFCluster& back = output.back();
0070     edm::Ref<std::vector<reco::PFRecHit> > seed;
0071     double energy = 0.0, highest_energy = 0.0;
0072     auto hitsAndFractions = sc.hits_and_fractions();
0073     for (const auto& hAndF : hitsAndFractions) {
0074       auto itr = detIdToIndex.find(hAndF.first);
0075       if (itr == detIdToIndex.end())
0076         continue;  // hit wasn't saved in reco
0077       auto ref = makeRefhit(input, itr->second);
0078       const double hit_energy = hAndF.second * ref->energy();
0079       energy += hit_energy;
0080       back.addRecHitFraction(reco::PFRecHitFraction(ref, hAndF.second));
0081       if (hit_energy > highest_energy || highest_energy == 0.0) {
0082         highest_energy = hit_energy;
0083         seed = ref;
0084       }
0085     }
0086     if (!back.hitsAndFractions().empty()) {
0087       back.setSeed(seed->detId());
0088       back.setEnergy(energy);
0089       back.setCorrectedEnergy(energy);
0090     } else {
0091       back.setSeed(-1);
0092       back.setEnergy(0.f);
0093     }
0094   }
0095 }