CMS 3D CMS Logo

HGCalClusteringAlgoBase.h
Go to the documentation of this file.
1 #ifndef RecoLocalCalo_HGCalRecProducers_HGCalClusteringAlgoBase_h
2 #define RecoLocalCalo_HGCalRecProducers_HGCalClusteringAlgoBase_h
3 
5 
9 
11 
12 // C/C++ headers
13 #include <vector>
14 #include <numeric>
15 
16 namespace hgcal_clustering {
17  template <typename T>
18  std::vector<size_t> sorted_indices(const std::vector<T> &v) {
19  // initialize original index locations
20  std::vector<size_t> idx(v.size());
21  std::iota(std::begin(idx), std::end(idx), 0);
22 
23  // sort indices based on comparing values in v
24  std::sort(idx.begin(), idx.end(), [&v](size_t i1, size_t i2) { return v[i1] > v[i2]; });
25 
26  return idx;
27  }
28 
29  template <typename T>
30  size_t max_index(const std::vector<T> &v) {
31  // initialize original index locations
32  std::vector<size_t> idx(v.size(), 0);
33  std::iota(std::begin(idx), std::end(idx), 0);
34 
35  // take the max index based on comparing values in v
36  auto maxidx = std::max_element(
37  idx.begin(), idx.end(), [&v](size_t i1, size_t i2) { return v[i1].data.rho < v[i2].data.rho; });
38 
39  return (*maxidx);
40  }
41 
42  //Density collection
43  typedef std::map<DetId, float> Density;
44 
45 }; // namespace hgcal_clustering
46 
48 public:
49  enum VerbosityLevel { pDEBUG = 0, pWARNING = 1, pINFO = 2, pERROR = 3 };
50 
53 
54  virtual void populate(const HGCRecHitCollection &hits) = 0;
55  virtual void makeClusters() = 0;
56  virtual std::vector<reco::BasicCluster> getClusters(bool) = 0;
57  virtual void reset() = 0;
58  virtual hgcal_clustering::Density getDensity() = 0;
59  virtual void getEventSetupPerAlgorithm(const edm::EventSetup &es) {}
60 
61  inline void getEventSetup(const edm::EventSetup &es) {
62  rhtools_.getEventSetup(es);
63  maxlayer_ = rhtools_.lastLayerBH();
64  lastLayerEE_ = rhtools_.lastLayerEE();
65  lastLayerFH_ = rhtools_.lastLayerFH();
66  firstLayerBH_ = rhtools_.firstLayerBH();
67  scintMaxIphi_ = rhtools_.getScintMaxIphi();
68  getEventSetupPerAlgorithm(es);
69  }
70  inline void setVerbosity(VerbosityLevel the_verbosity) { verbosity_ = the_verbosity; }
71  inline void setAlgoId(reco::CaloCluster::AlgoId algo) { algoId_ = algo; }
72 
73  //max number of layers
74  unsigned int maxlayer_;
75  // last layer per subdetector
76  unsigned int lastLayerEE_;
77  unsigned int lastLayerFH_;
78  unsigned int firstLayerBH_;
80 
81 protected:
82  // The verbosity level
84 
85  // The vector of clusters
86  std::vector<reco::BasicCluster> clusters_v_;
87 
89 
90  // The algo id
92 };
93 
94 #endif
std::vector< reco::BasicCluster > clusters_v_
void setVerbosity(VerbosityLevel the_verbosity)
std::map< DetId, float > Density
HGCalClusteringAlgoBase(VerbosityLevel v, reco::CaloCluster::AlgoId algo)
std::vector< size_t > sorted_indices(const std::vector< T > &v)
#define end
Definition: vmac.h:39
reco::CaloCluster::AlgoId algoId_
#define begin
Definition: vmac.h:32
void getEventSetup(const edm::EventSetup &es)
size_t max_index(const std::vector< T > &v)
void setAlgoId(reco::CaloCluster::AlgoId algo)
void reset(double vett[256])
Definition: TPedValues.cc:11
virtual void getEventSetupPerAlgorithm(const edm::EventSetup &es)