Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-06-22 02:24:04

0001 #include "FWCore/Framework/interface/ConsumesCollector.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/Frameworkfwd.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/Framework/interface/stream/EDProducer.h"
0006 #include "RecoParticleFlow/PFClusterProducer/interface/InitialClusteringStepBase.h"
0007 #include "RecoParticleFlow/PFClusterProducer/interface/PFCPositionCalculatorBase.h"
0008 #include "RecoParticleFlow/PFClusterProducer/interface/PFClusterBuilderBase.h"
0009 #include "RecoParticleFlow/PFClusterProducer/interface/PFClusterEnergyCorrectorBase.h"
0010 #include "RecoParticleFlow/PFClusterProducer/interface/RecHitTopologicalCleanerBase.h"
0011 #include "RecoParticleFlow/PFClusterProducer/interface/SeedFinderBase.h"
0012 #include "CondFormats/DataRecord/interface/HcalPFCutsRcd.h"
0013 #include "CondTools/Hcal/interface/HcalPFCutsHandler.h"
0014 
0015 #include <memory>
0016 
0017 class PFMultiDepthClusterProducer : public edm::stream::EDProducer<> {
0018   typedef InitialClusteringStepBase ICSB;
0019   typedef PFClusterBuilderBase PFCBB;
0020   typedef PFCPositionCalculatorBase PosCalc;
0021 
0022 public:
0023   PFMultiDepthClusterProducer(const edm::ParameterSet&);
0024   ~PFMultiDepthClusterProducer() override = default;
0025 
0026   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0027   void produce(edm::Event&, const edm::EventSetup&) override;
0028   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0029 
0030 private:
0031   // inputs
0032   edm::EDGetTokenT<reco::PFClusterCollection> _clustersLabel;
0033   edm::ESGetToken<HcalPFCuts, HcalPFCutsRcd> hcalCutsToken_;
0034   // options
0035   // the actual algorithm
0036   std::unique_ptr<PFClusterBuilderBase> _pfClusterBuilder;
0037   std::unique_ptr<PFClusterEnergyCorrectorBase> _energyCorrector;
0038   bool cutsFromDB;
0039   HcalPFCuts const* paramPF = nullptr;
0040 };
0041 
0042 DEFINE_FWK_MODULE(PFMultiDepthClusterProducer);
0043 
0044 void PFMultiDepthClusterProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0045   edm::ParameterSetDescription desc;
0046   desc.add<edm::InputTag>("clustersSource", {});
0047   desc.add<edm::ParameterSetDescription>("energyCorrector", {});
0048   {
0049     edm::ParameterSetDescription pset0;
0050     pset0.add<std::string>("algoName", "PFMultiDepthClusterizer");
0051     {
0052       edm::ParameterSetDescription pset1;
0053       pset1.add<std::string>("algoName", "Basic2DGenericPFlowPositionCalc");
0054       {
0055         edm::ParameterSetDescription psd;
0056         psd.add<std::vector<int>>("depths", {});
0057         psd.add<std::string>("detector", "");
0058         psd.add<std::vector<double>>("logWeightDenominator", {});
0059         pset1.addVPSet("logWeightDenominatorByDetector", psd, {});
0060       }
0061       pset1.add<double>("minAllowedNormalization", 1e-09);
0062       pset1.add<double>("minFractionInCalc", 1e-09);
0063       pset1.add<int>("posCalcNCrystals", -1);
0064       pset1.add<edm::ParameterSetDescription>("timeResolutionCalcBarrel", {});
0065       pset1.add<edm::ParameterSetDescription>("timeResolutionCalcEndcap", {});
0066       pset0.add<edm::ParameterSetDescription>("allCellsPositionCalc", pset1);
0067     }
0068     pset0.add<edm::ParameterSetDescription>("positionCalc", {});
0069     pset0.add<double>("minFractionToKeep", 1e-07);
0070     pset0.add<double>("nSigmaEta", 2.0);
0071     pset0.add<double>("nSigmaPhi", 2.0);
0072     desc.add<edm::ParameterSetDescription>("pfClusterBuilder", pset0);
0073   }
0074   desc.add<edm::ParameterSetDescription>("positionReCalc", {});
0075   desc.add<bool>("usePFThresholdsFromDB", false);
0076   descriptions.addWithDefaultLabel(desc);
0077 }
0078 
0079 #ifdef PFLOW_DEBUG
0080 #define LOGVERB(x) edm::LogVerbatim(x)
0081 #define LOGWARN(x) edm::LogWarning(x)
0082 #define LOGERR(x) edm::LogError(x)
0083 #define LOGDRESSED(x) edm::LogInfo(x)
0084 #else
0085 #define LOGVERB(x) LogTrace(x)
0086 #define LOGWARN(x) edm::LogWarning(x)
0087 #define LOGERR(x) edm::LogError(x)
0088 #define LOGDRESSED(x) LogDebug(x)
0089 #endif
0090 
0091 PFMultiDepthClusterProducer::PFMultiDepthClusterProducer(const edm::ParameterSet& conf) {
0092   _clustersLabel = consumes<reco::PFClusterCollection>(conf.getParameter<edm::InputTag>("clustersSource"));
0093   cutsFromDB = conf.getParameter<bool>("usePFThresholdsFromDB");
0094   const edm::ParameterSet& pfcConf = conf.getParameterSet("pfClusterBuilder");
0095 
0096   edm::ConsumesCollector&& cc = consumesCollector();
0097 
0098   if (cutsFromDB) {
0099     hcalCutsToken_ = esConsumes<HcalPFCuts, HcalPFCutsRcd, edm::Transition::BeginRun>(edm::ESInputTag("", "withTopo"));
0100   }
0101 
0102   if (!pfcConf.empty()) {
0103     const std::string& pfcName = pfcConf.getParameter<std::string>("algoName");
0104     if (!pfcName.empty())
0105       _pfClusterBuilder = PFClusterBuilderFactory::get()->create(pfcName, pfcConf, cc);
0106   }
0107   // see if new need to apply corrections, setup if there.
0108   const edm::ParameterSet& cConf = conf.getParameterSet("energyCorrector");
0109   if (!cConf.empty()) {
0110     const std::string& cName = cConf.getParameter<std::string>("algoName");
0111     if (!cName.empty())
0112       _energyCorrector = PFClusterEnergyCorrectorFactory::get()->create(cName, cConf);
0113   }
0114 
0115   produces<reco::PFClusterCollection>();
0116 }
0117 
0118 void PFMultiDepthClusterProducer::beginRun(const edm::Run& run, const edm::EventSetup& es) {
0119   if (cutsFromDB) {
0120     paramPF = &es.getData(hcalCutsToken_);
0121   }
0122   _pfClusterBuilder->update(es);
0123 }
0124 
0125 void PFMultiDepthClusterProducer::produce(edm::Event& e, const edm::EventSetup& es) {
0126   _pfClusterBuilder->reset();
0127 
0128   edm::Handle<reco::PFClusterCollection> inputClusters;
0129   e.getByToken(_clustersLabel, inputClusters);
0130 
0131   std::vector<bool> seedable;
0132 
0133   auto pfClusters = std::make_unique<reco::PFClusterCollection>();
0134   _pfClusterBuilder->buildClusters(*inputClusters, seedable, *pfClusters, paramPF);
0135   LOGVERB("PFMultiDepthClusterProducer::produce()") << *_pfClusterBuilder;
0136 
0137   if (_energyCorrector) {
0138     _energyCorrector->correctEnergies(*pfClusters);
0139   }
0140   e.put(std::move(pfClusters));
0141 }