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 
16 namespace edm {
17  class Event;
18  class EventSetup;
19 }
20 
22 
25  public:
27  edm::ConsumesCollector& sumes):
28  _nSeeds(0), _nClustersFound(0),
29  _layerMap({ {"PS2",(int)PFLayer::PS2},
30  {"PS1",(int)PFLayer::PS1},
31  {"ECAL_ENDCAP",(int)PFLayer::ECAL_ENDCAP},
32  {"ECAL_BARREL",(int)PFLayer::ECAL_BARREL},
33  {"NONE",(int)PFLayer::NONE},
34  {"HCAL_BARREL1",(int)PFLayer::HCAL_BARREL1},
35  {"HCAL_BARREL2_RING0",(int)PFLayer::HCAL_BARREL2},
36  {"HCAL_BARREL2_RING1",100*(int)PFLayer::HCAL_BARREL2},
37  {"HCAL_ENDCAP",(int)PFLayer::HCAL_ENDCAP},
38  {"HF_EM",(int)PFLayer::HF_EM},
39  {"HF_HAD",(int)PFLayer::HF_HAD},
40  {"HGCAL",(int)PFLayer::HGCAL} }),
41  _algoName(conf.getParameter<std::string>("algoName")) {
42  const std::vector<edm::ParameterSet>& thresholds =
43  conf.getParameterSetVector("thresholdsByDetector");
44  for( const auto& pset : thresholds ) {
45  const std::string& det = pset.getParameter<std::string>("detector");
46  const double& thresh_E =
47  pset.getParameter<double>("gatheringThreshold");
48  const double& thresh_pT =
49  pset.getParameter<double>("gatheringThresholdPt");
50  const double thresh_pT2 = thresh_pT*thresh_pT;
51  auto entry = _layerMap.find(det);
52  if( entry == _layerMap.end() ) {
53  throw cms::Exception("InvalidDetectorLayer")
54  << "Detector layer : " << det << " is not in the list of recognized"
55  << " detector layers!";
56  }
57  _thresholds.emplace(_layerMap.find(det)->second,
58  std::make_pair(thresh_E,thresh_pT2));
59  }
60  }
61  virtual ~InitialClusteringStepBase() = default;
62  // get rid of things we should never use...
63  InitialClusteringStepBase(const ICSB&) = delete;
64  ICSB& operator=(const ICSB&) = delete;
65 
66  virtual void update(const edm::EventSetup&) { }
67 
68  virtual void updateEvent(const edm::Event&) { }
69 
70  virtual void buildClusters(const edm::Handle<reco::PFRecHitCollection>&,
71  const std::vector<bool>& mask, // mask flags
72  const std::vector<bool>& seeds, // seed flags
73  reco::PFClusterCollection&) = 0; //output
74 
75  std::ostream& operator<<(std::ostream& o) {
76  o << "InitialClusteringStep with algo \"" << _algoName
77  << "\" located " << _nSeeds << " seeds and built "
78  << _nClustersFound << " clusters from those seeds. ";
79  return o;
80  }
81 
82  void reset() { _nSeeds = _nClustersFound = 0; }
83 
84  protected:
86  const unsigned i ) const {
87  return reco::PFRecHitRef(h,i);
88  }
89  unsigned _nSeeds, _nClustersFound; // basic performance information
90  const std::unordered_map<std::string,int> _layerMap;
91  std::unordered_map<int,std::pair<double,double> >
93 
94  private:
96 
97 };
98 
101 
102 #endif
const std::unordered_map< std::string, int > _layerMap
std::ostream & operator<<(std::ostream &o)
reco::PFRecHitRef makeRefhit(const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
InitialClusteringStepBase ICSB
edmplugin::PluginFactory< InitialClusteringStepBase *(const edm::ParameterSet &, edm::ConsumesCollector &) > InitialClusteringStepFactory
virtual void updateEvent(const edm::Event &)
edm::Ref< PFRecHitCollection > PFRecHitRef
persistent reference to PFRecHit objects
Definition: PFRecHitFwd.h:15
virtual void update(const edm::EventSetup &)
InitialClusteringStepBase(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
HLT enums.
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
std::unordered_map< int, std::pair< double, double > > _thresholds