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 }
21 
23 
26  public:
28  edm::ConsumesCollector& sumes):
29  _nSeeds(0), _nClustersFound(0),
30  _layerMap({ {"PS2",(int)PFLayer::PS2},
31  {"PS1",(int)PFLayer::PS1},
32  {"ECAL_ENDCAP",(int)PFLayer::ECAL_ENDCAP},
33  {"ECAL_BARREL",(int)PFLayer::ECAL_BARREL},
34  {"NONE",(int)PFLayer::NONE},
35  {"HCAL_BARREL1",(int)PFLayer::HCAL_BARREL1},
36  {"HCAL_BARREL2_RING0",(int)PFLayer::HCAL_BARREL2},
37  {"HCAL_BARREL2_RING1",100*(int)PFLayer::HCAL_BARREL2},
38  {"HCAL_ENDCAP",(int)PFLayer::HCAL_ENDCAP},
39  {"HF_EM",(int)PFLayer::HF_EM},
40  {"HF_HAD",(int)PFLayer::HF_HAD},
41  {"HGCAL",(int)PFLayer::HGCAL} }),
42  _algoName(conf.getParameter<std::string>("algoName")) {
43  const std::vector<edm::ParameterSet>& thresholds =
44  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,
78  std::make_tuple(depths,thresh_E,thresh_pT2));
79  }
80  }
81  virtual ~InitialClusteringStepBase() = default;
82  // get rid of things we should never use...
83  InitialClusteringStepBase(const ICSB&) = delete;
84  ICSB& operator=(const ICSB&) = delete;
85 
86  virtual void update(const edm::EventSetup&) { }
87 
88  virtual void updateEvent(const edm::Event&) { }
89 
90  virtual void buildClusters(const edm::Handle<reco::PFRecHitCollection>&,
91  const std::vector<bool>& mask, // mask flags
92  const std::vector<bool>& seeds, // seed flags
93  reco::PFClusterCollection&) = 0; //output
94 
95  std::ostream& operator<<(std::ostream& o) const {
96  o << "InitialClusteringStep with algo \"" << _algoName
97  << "\" located " << _nSeeds << " seeds and built "
98  << _nClustersFound << " clusters from those seeds. ";
99  return o;
100  }
101 
102  void reset() { _nSeeds = _nClustersFound = 0; }
103 
104  protected:
106  const unsigned i ) const {
107  return reco::PFRecHitRef(h,i);
108  }
109  unsigned _nSeeds, _nClustersFound; // basic performance information
110  const std::unordered_map<std::string,int> _layerMap;
111 
112  typedef std::tuple<std::vector<int> ,std::vector<double> , std::vector<double> > I3tuple;
113  std::unordered_map<int,I3tuple > _thresholds;
114 
115  private:
117 
118 };
119 
120 std::ostream& operator<<(std::ostream& o, const InitialClusteringStepBase& a);
121 
124 
125 #endif
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
std::unordered_map< int, I3tuple > _thresholds
const std::unordered_map< std::string, int > _layerMap
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 &)
std::tuple< std::vector< int >,std::vector< double >, std::vector< double > > I3tuple
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.
double a
Definition: hdecay.h:121
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
std::ostream & operator<<(std::ostream &ost, const HLTGlobalStatus &hlt)
Formatted printout of trigger tbale.
std::ostream & operator<<(std::ostream &o) const