File indexing completed on 2023-03-17 11:21:00
0001 #ifndef __InitialClusteringStepBase_H__
0002 #define __InitialClusteringStepBase_H__
0003
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0007 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFwd.h"
0008 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0009 #include "DataFormats/ParticleFlowReco/interface/PFClusterFwd.h"
0010 #include "DataFormats/Common/interface/RefToBaseVector.h"
0011
0012 #include <string>
0013 #include <iostream>
0014 #include <unordered_map>
0015 #include <tuple>
0016
0017 namespace edm {
0018 class Event;
0019 class EventSetup;
0020 }
0021
0022 #include "FWCore/Framework/interface/ConsumesCollector.h"
0023
0024 class InitialClusteringStepBase {
0025 typedef InitialClusteringStepBase ICSB;
0026
0027 public:
0028 InitialClusteringStepBase(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0029 : _nSeeds(0),
0030 _nClustersFound(0),
0031 _layerMap({{"PS2", (int)PFLayer::PS2},
0032 {"PS1", (int)PFLayer::PS1},
0033 {"ECAL_ENDCAP", (int)PFLayer::ECAL_ENDCAP},
0034 {"ECAL_BARREL", (int)PFLayer::ECAL_BARREL},
0035 {"NONE", (int)PFLayer::NONE},
0036 {"HCAL_BARREL1", (int)PFLayer::HCAL_BARREL1},
0037 {"HCAL_BARREL2_RING0", (int)PFLayer::HCAL_BARREL2},
0038 {"HCAL_BARREL2_RING1", 100 * (int)PFLayer::HCAL_BARREL2},
0039 {"HCAL_ENDCAP", (int)PFLayer::HCAL_ENDCAP},
0040 {"HF_EM", (int)PFLayer::HF_EM},
0041 {"HF_HAD", (int)PFLayer::HF_HAD},
0042 {"HGCAL", (int)PFLayer::HGCAL}}),
0043 _algoName(conf.getParameter<std::string>("algoName")) {
0044 const std::vector<edm::ParameterSet>& thresholds = conf.getParameterSetVector("thresholdsByDetector");
0045 for (const auto& pset : thresholds) {
0046 const std::string& det = pset.getParameter<std::string>("detector");
0047
0048 std::vector<int> depths;
0049 std::vector<double> thresh_E;
0050 std::vector<double> thresh_pT;
0051 std::vector<double> thresh_pT2;
0052
0053 if (det == std::string("HCAL_BARREL1") || det == std::string("HCAL_ENDCAP")) {
0054 depths = pset.getParameter<std::vector<int> >("depths");
0055 thresh_E = pset.getParameter<std::vector<double> >("gatheringThreshold");
0056 thresh_pT = pset.getParameter<std::vector<double> >("gatheringThresholdPt");
0057 if (thresh_E.size() != depths.size() || thresh_pT.size() != depths.size()) {
0058 throw cms::Exception("InvalidGatheringThreshold")
0059 << "gatheringThresholds mismatch with the numbers of depths";
0060 }
0061 } else {
0062 depths.push_back(0);
0063 thresh_E.push_back(pset.getParameter<double>("gatheringThreshold"));
0064 thresh_pT.push_back(pset.getParameter<double>("gatheringThresholdPt"));
0065 }
0066
0067 for (unsigned int i = 0; i < thresh_pT.size(); ++i) {
0068 thresh_pT2.push_back(thresh_pT[i] * thresh_pT[i]);
0069 }
0070
0071 auto entry = _layerMap.find(det);
0072 if (entry == _layerMap.end()) {
0073 throw cms::Exception("InvalidDetectorLayer")
0074 << "Detector layer : " << det << " is not in the list of recognized"
0075 << " detector layers!";
0076 }
0077 _thresholds.emplace(_layerMap.find(det)->second, std::make_tuple(depths, thresh_E, thresh_pT2));
0078 }
0079 }
0080 virtual ~InitialClusteringStepBase() = default;
0081
0082 InitialClusteringStepBase(const ICSB&) = delete;
0083 ICSB& operator=(const ICSB&) = delete;
0084
0085 virtual void update(const edm::EventSetup&) {}
0086
0087 virtual void updateEvent(const edm::Event&) {}
0088
0089 virtual void buildClusters(const edm::Handle<reco::PFRecHitCollection>&,
0090 const std::vector<bool>& mask,
0091 const std::vector<bool>& seeds,
0092 reco::PFClusterCollection&) = 0;
0093
0094 std::ostream& operator<<(std::ostream& o) const {
0095 o << "InitialClusteringStep with algo \"" << _algoName << "\" located " << _nSeeds << " seeds and built "
0096 << _nClustersFound << " clusters from those seeds. ";
0097 return o;
0098 }
0099
0100 void reset() { _nSeeds = _nClustersFound = 0; }
0101
0102 protected:
0103 reco::PFRecHitRef makeRefhit(const edm::Handle<reco::PFRecHitCollection>& h, const unsigned i) const {
0104 return reco::PFRecHitRef(h, i);
0105 }
0106 unsigned _nSeeds, _nClustersFound;
0107 const std::unordered_map<std::string, int> _layerMap;
0108
0109 typedef std::tuple<std::vector<int>, std::vector<double>, std::vector<double> > I3tuple;
0110 std::unordered_map<int, I3tuple> _thresholds;
0111
0112 private:
0113 const std::string _algoName;
0114 };
0115
0116 std::ostream& operator<<(std::ostream& o, const InitialClusteringStepBase& a);
0117
0118 #include "FWCore/PluginManager/interface/PluginFactory.h"
0119 typedef edmplugin::PluginFactory<InitialClusteringStepBase*(const edm::ParameterSet&, edm::ConsumesCollector&)>
0120 InitialClusteringStepFactory;
0121
0122 #endif