CMS 3D CMS Logo

InitialClusteringStepBase.h
Go to the documentation of this file.
1 #ifndef __InitialClusteringStepBase_H__
2 #define __InitialClusteringStepBase_H__
3 
11 
12 #include <string>
13 #include <iostream>
14 #include <unordered_map>
15 #include <tuple>
16 
17 namespace edm {
18  class Event;
19  class EventSetup;
20 } // namespace edm
21 
23 
26 
27 public:
29  : _nSeeds(0),
30  _nClustersFound(0),
31  _layerMap({{"PS2", (int)PFLayer::PS2},
32  {"PS1", (int)PFLayer::PS1},
33  {"ECAL_ENDCAP", (int)PFLayer::ECAL_ENDCAP},
34  {"ECAL_BARREL", (int)PFLayer::ECAL_BARREL},
35  {"NONE", (int)PFLayer::NONE},
36  {"HCAL_BARREL1", (int)PFLayer::HCAL_BARREL1},
37  {"HCAL_BARREL2_RING0", (int)PFLayer::HCAL_BARREL2},
38  {"HCAL_BARREL2_RING1", 100 * (int)PFLayer::HCAL_BARREL2},
39  {"HCAL_ENDCAP", (int)PFLayer::HCAL_ENDCAP},
40  {"HF_EM", (int)PFLayer::HF_EM},
41  {"HF_HAD", (int)PFLayer::HF_HAD},
42  {"HGCAL", (int)PFLayer::HGCAL}}),
43  _algoName(conf.getParameter<std::string>("algoName")) {
44  const std::vector<edm::ParameterSet>& thresholds = conf.getParameterSetVector("thresholdsByDetector");
45  for (const auto& pset : thresholds) {
46  const std::string& det = pset.getParameter<std::string>("detector");
47 
48  std::vector<int> depths;
49  std::vector<double> thresh_E;
50  std::vector<double> thresh_pT;
51  std::vector<double> thresh_pT2;
52 
53  if (det == std::string("HCAL_BARREL1") || det == std::string("HCAL_ENDCAP")) {
54  depths = pset.getParameter<std::vector<int> >("depths");
55  thresh_E = pset.getParameter<std::vector<double> >("gatheringThreshold");
56  thresh_pT = pset.getParameter<std::vector<double> >("gatheringThresholdPt");
57  if (thresh_E.size() != depths.size() || thresh_pT.size() != depths.size()) {
58  throw cms::Exception("InvalidGatheringThreshold")
59  << "gatheringThresholds mismatch with the numbers of depths";
60  }
61  } else {
62  depths.push_back(0);
63  thresh_E.push_back(pset.getParameter<double>("gatheringThreshold"));
64  thresh_pT.push_back(pset.getParameter<double>("gatheringThresholdPt"));
65  }
66 
67  for (unsigned int i = 0; i < thresh_pT.size(); ++i) {
68  thresh_pT2.push_back(thresh_pT[i] * thresh_pT[i]);
69  }
70 
71  auto entry = _layerMap.find(det);
72  if (entry == _layerMap.end()) {
73  throw cms::Exception("InvalidDetectorLayer")
74  << "Detector layer : " << det << " is not in the list of recognized"
75  << " detector layers!";
76  }
77  _thresholds.emplace(_layerMap.find(det)->second, std::make_tuple(depths, thresh_E, thresh_pT2));
78  }
79  }
80  virtual ~InitialClusteringStepBase() = default;
81  // get rid of things we should never use...
82  InitialClusteringStepBase(const ICSB&) = delete;
83  ICSB& operator=(const ICSB&) = delete;
84 
85  virtual void update(const edm::EventSetup&) {}
86 
87  virtual void updateEvent(const edm::Event&) {}
88 
90  const std::vector<bool>& mask, // mask flags
91  const std::vector<bool>& seeds, // seed flags
92  reco::PFClusterCollection&) = 0; //output
93 
94  std::ostream& operator<<(std::ostream& o) const {
95  o << "InitialClusteringStep with algo \"" << _algoName << "\" located " << _nSeeds << " seeds and built "
96  << _nClustersFound << " clusters from those seeds. ";
97  return o;
98  }
99 
100  void reset() { _nSeeds = _nClustersFound = 0; }
101 
102 protected:
104  return reco::PFRecHitRef(h, i);
105  }
106  unsigned _nSeeds, _nClustersFound; // basic performance information
107  const std::unordered_map<std::string, int> _layerMap;
108 
109  typedef std::tuple<std::vector<int>, std::vector<double>, std::vector<double> > I3tuple;
110  std::unordered_map<int, I3tuple> _thresholds;
111 
112 private:
114 };
115 
116 std::ostream& operator<<(std::ostream& o, const InitialClusteringStepBase& a);
117 
121 
122 #endif
const std::unordered_map< std::string, int > _layerMap
reco::PFRecHitRef makeRefhit(const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
std::ostream & operator<<(std::ostream &o) const
constexpr uint32_t mask
Definition: gpuClustering.h:26
InitialClusteringStepBase ICSB
virtual void updateEvent(const edm::Event &)
virtual void buildClusters(const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &mask, const std::vector< bool > &seeds, reco::PFClusterCollection &)=0
InitialClusteringStepBase(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
edm::Ref< PFRecHitCollection > PFRecHitRef
persistent reference to PFRecHit objects
Definition: PFRecHitFwd.h:15
edmplugin::PluginFactory< InitialClusteringStepBase *(const edm::ParameterSet &, edm::ConsumesCollector &)> InitialClusteringStepFactory
virtual void update(const edm::EventSetup &)
std::tuple< std::vector< int >, std::vector< double >, std::vector< double > > I3tuple
std::ostream & operator<<(std::ostream &o, const InitialClusteringStepBase &a)
ICSB & operator=(const ICSB &)=delete
HLT enums.
double a
Definition: hdecay.h:121
virtual ~InitialClusteringStepBase()=default
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
std::unordered_map< int, I3tuple > _thresholds