CMS 3D CMS Logo

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