CMS 3D CMS Logo

CaloClusterer.h
Go to the documentation of this file.
1 #ifndef L1Trigger_Phase2L1ParticleFlow_CALOCLUSTERER_H
2 #define L1Trigger_Phase2L1ParticleFlow_CALOCLUSTERER_H
3 
7 // fwd declarations
8 namespace edm {
9  class ParameterSet;
10 }
11 
12 // real includes
13 #include <cstdint>
14 #include <cmath>
15 #include <vector>
16 #include <array>
17 #include <algorithm>
20 
21 namespace l1tpf_calo {
22  class Grid {
23  public:
24  virtual ~Grid() {}
25  unsigned int size() const { return ncells_; }
26  virtual int find_cell(float eta, float phi) const = 0;
27  int neighbour(int icell, unsigned int idx) const { return neighbours_[icell][idx]; }
28  float eta(int icell) const { return eta_[icell]; }
29  float phi(int icell) const { return phi_[icell]; }
30  float etaWidth(int icell) const { return etaWidth_[icell]; }
31  float phiWidth(int icell) const { return phiWidth_[icell]; }
32  int ieta(int icell) const { return ieta_[icell]; }
33  int iphi(int icell) const { return iphi_[icell]; }
34 
35  protected:
36  Grid(unsigned int size)
37  : ncells_(size),
38  eta_(size),
39  etaWidth_(size),
40  phi_(size),
41  phiWidth_(size),
42  ieta_(size),
43  iphi_(size),
44  neighbours_(size) {}
45  unsigned int ncells_;
46  std::vector<float> eta_, etaWidth_, phi_, phiWidth_;
47  std::vector<int> ieta_, iphi_;
48  std::vector<std::array<int, 8>> neighbours_; // indices of the neigbours, -1 = none
49  };
50 
51  class Phase1GridBase : public Grid {
52  public:
53  Phase1GridBase(int nEta, int nPhi, int ietaCoarse, int ietaVeryCoarse, const float *towerEtas);
54 
55  int find_cell(float eta, float phi) const override;
56  int ifind_cell(int ieta, int iphi) const { return cell_map_[(ieta + nEta_) + 2 * nEta_ * (iphi - 1)]; }
57 
58  protected:
60  const float *towerEtas_;
61  std::vector<int> cell_map_;
62  // valid ieta, iphi (does not check for outside bounds, only for non-existence of ieta=0, iphi=0, and coarser towers at high eta)
63  bool valid_ieta_iphi(int ieta, int iphi) const {
64  if (ieta == 0 || iphi == 0)
65  return false;
66  if (std::abs(ieta) >= ietaVeryCoarse_ && (iphi % 4 != 1))
67  return false;
68  if (std::abs(ieta) >= ietaCoarse_ && (iphi % 2 != 1))
69  return false;
70  return true;
71  }
72  // move by +/-1 around a cell; return icell or -1 if not available
73  int imove(int ieta, int iphi, int deta, int dphi);
74  };
75 
76  class Phase1Grid : public Phase1GridBase {
77  public:
80 
81  protected:
82  static const int phase1_nEta_ = 41, phase1_nPhi_ = 72, phase1_ietaCoarse_ = 29, phase1_ietaVeryCoarse_ = 40;
83  static const float phase1_towerEtas_[phase1_nEta_];
84  };
85  class Phase2Grid : public Phase1GridBase {
86  public:
89 
90  protected:
91  static const int phase2_nEta_ = 48, phase2_nPhi_ = 72, phase2_ietaCoarse_ = 36, phase2_ietaVeryCoarse_ = 47;
92  static const float phase2_towerEtas_[phase2_nEta_];
93  };
94 
95  template <typename T>
96  class GridData {
97  public:
98  GridData() : grid_(nullptr), data_(), empty_() {}
99  GridData(const Grid &grid) : grid_(&grid), data_(grid.size()), empty_() {}
100 
101  T &operator()(float eta, float phi) { return data_[grid_->find_cell(eta, phi)]; }
102  const T &operator()(float eta, float phi) const { return data_[grid_->find_cell(eta, phi)]; }
103 
104  const Grid &grid() const { return *grid_; }
105 
106  unsigned int size() const { return data_.size(); }
107 
108  float eta(int icell) const { return grid().eta(icell); }
109  float phi(int icell) const { return grid().phi(icell); }
110  int ieta(int icell) const { return grid().ieta(icell); }
111  int iphi(int icell) const { return grid().iphi(icell); }
112 
113  T &operator[](int icell) { return data_[icell]; }
114  const T &operator[](int icell) const { return data_[icell]; }
115 
116  const T &neigh(int icell, unsigned int idx) const {
117  int ineigh = grid_->neighbour(icell, idx);
118  return (ineigh < 0 ? empty_ : data_[ineigh]);
119  }
120 
122  assert(grid_ == other.grid_);
123  data_ = other.data_;
124  return *this;
125  }
127  assert(grid_ == other.grid_);
128  for (unsigned int i = 0, n = data_.size(); i < n; ++i) {
129  data_[i] += other.data_[i];
130  }
131  return *this;
132  }
133 
134  // always defined
135  void fill(const T &val) { std::fill(data_.begin(), data_.end(), val); }
136  void zero() { fill(T()); }
137 
138  // defined only if T has a 'clear' method
139  void clear() {
140  for (T &t : data_)
141  t.clear();
142  }
143 
144  private:
145  const Grid *grid_;
146  std::vector<T> data_;
147  const T empty_;
148  };
151 
152  struct PreCluster {
154  float ptLocalMax; // pt if it's a local max, zero otherwise
155  float ptOverNeighLocalMaxSum; // pt / (sum of ptLocalMax of neighbours); zero if no neighbours
157  };
159 
160  struct Cluster {
161  Cluster() : et(0), eta(0), phi(0) {}
162  float et, eta, phi;
163  std::vector<std::pair<int, float>> constituents;
164  void clear() {
165  et = eta = phi = 0;
166  constituents.clear();
167  }
168  };
169 
170  struct CombinedCluster : public Cluster {
171  float ecal_et, hcal_et;
172  void clear() {
173  Cluster::clear();
174  ecal_et = hcal_et = 0;
175  }
176  };
177 
178  const Grid *getGrid(const std::string &type);
179 
181  public:
184  void clear();
185  void add(const reco::Candidate &c) { add(c.pt(), c.eta(), c.phi()); }
186  void add(float pt, float eta, float phi) { rawet_(eta, phi) += pt; }
187  void run();
188 
190  // note: there can be some double-counting as the same unclustered energy can go into more clusters
191  void grow();
192 
193  const EtGrid &raw() const { return rawet_; }
194  const IndexGrid &indexGrid() const { return clusterIndex_; }
195  const std::vector<Cluster> &clusters() const { return clusters_; }
196  const Cluster &cluster(int i) const {
197  return (i == -1 || clusterIndex_[i] == -1) ? nullCluster_ : clusters_[clusterIndex_[i]];
198  }
199 
201  EtGrid &raw() { return rawet_; }
202 
203  // for the moment, generic interface that takes a cluster and returns the corrected pt
204  template <typename Corrector>
205  void correct(const Corrector &corrector) {
206  for (Cluster &c : clusters_) {
207  c.et = corrector(c);
208  }
209  }
210 
211  std::unique_ptr<l1t::PFClusterCollection> fetchCells(bool unclusteredOnly = false, float ptMin = 0.) const;
212 
213  std::unique_ptr<l1t::PFClusterCollection> fetch(float ptMin = 0.) const;
214  std::unique_ptr<l1t::PFClusterCollection> fetch(const edm::OrphanHandle<l1t::PFClusterCollection> &cells,
215  float ptMin = 0.) const;
216 
217  private:
218  enum class EnergyShareAlgo {
219  Fractions, /* each local maximum neighbour takes a share proportional to its value */
220  None, /* each local maximum neighbour takes all the value (double counting!) */
221  Greedy, /* assing cell to the highest local maximum neighbour */
222  Crude
223  }; /* if there's more than one local maximum neighbour, they all take half of the value (no fp division) */
224  const Grid *grid_;
228  std::vector<Cluster> clusters_;
232  bool energyWeightedPosition_; // do the energy-weighted cluster position instead of the cell center
233  };
234 
236  public:
238  const SingleCaloClusterer &ecal,
239  const SingleCaloClusterer &hcal);
240  virtual ~SimpleCaloLinkerBase();
241  virtual void clear() { clearBase(); }
242  virtual void run() = 0;
243  void clearBase() {
244  clusters_.clear();
245  clusterIndex_.fill(-1);
246  }
247 
248  // for the moment, generic interface that takes a cluster and returns the corrected pt
249  template <typename Corrector>
250  void correct(const Corrector &corrector) {
251  for (CombinedCluster &c : clusters_) {
252  c.et = corrector(c);
253  }
254  }
255 
256  std::unique_ptr<l1t::PFClusterCollection> fetch() const;
257  std::unique_ptr<l1t::PFClusterCollection> fetch(const edm::OrphanHandle<l1t::PFClusterCollection> &ecal,
259 
260  protected:
261  const Grid *grid_;
264  std::vector<CombinedCluster> clusters_;
267  };
268 
270  public:
272  ~SimpleCaloLinker() override;
273  void clear() override;
274  void run() override;
275 
276  protected:
278  };
280  public:
282  ~FlatCaloLinker() override;
283  void clear() override;
284  void run() override;
285 
286  protected:
288  };
289 
290  // makes a calo linker (pointer will be owned by the callee)
291  std::unique_ptr<SimpleCaloLinkerBase> makeCaloLinker(const edm::ParameterSet &pset,
292  const SingleCaloClusterer &ecal,
293  const SingleCaloClusterer &hcal);
294 
295 } // namespace l1tpf_calo
296 
297 #endif
Phase1GridBase(int nEta, int nPhi, int ietaCoarse, int ietaVeryCoarse, const float *towerEtas)
std::vector< CombinedCluster > clusters_
int iphi(int icell) const
Definition: CaloClusterer.h:33
GridData< float > EtGrid
int imove(int ieta, int iphi, int deta, int dphi)
float eta(int icell) const
Definition: CaloClusterer.h:28
std::vector< float > eta_
Definition: CaloClusterer.h:46
const Grid & grid() const
T & operator()(float eta, float phi)
int ieta(int icell) const
Definition: CaloClusterer.h:32
void fill(const T &val)
SimpleCaloLinker(const edm::ParameterSet &pset, const SingleCaloClusterer &ecal, const SingleCaloClusterer &hcal)
const SingleCaloClusterer & hcal_
SingleCaloClusterer(const edm::ParameterSet &pset)
constexpr float ptMin
const IndexGrid & indexGrid() const
int iphi(int icell) const
std::unique_ptr< SimpleCaloLinkerBase > makeCaloLinker(const edm::ParameterSet &pset, const SingleCaloClusterer &ecal, const SingleCaloClusterer &hcal)
assert(be >=bs)
int ieta(int icell) const
GridData< PreCluster > PreClusterGrid
void add(const reco::Candidate &c)
unsigned int ncells_
Definition: CaloClusterer.h:45
static const int phase1_ietaVeryCoarse_
Definition: CaloClusterer.h:82
std::vector< float > etaWidth_
Definition: CaloClusterer.h:46
std::vector< Cluster > clusters_
const T & operator[](int icell) const
FlatCaloLinker(const edm::ParameterSet &pset, const SingleCaloClusterer &ecal, const SingleCaloClusterer &hcal)
SimpleCaloLinkerBase(const edm::ParameterSet &pset, const SingleCaloClusterer &ecal, const SingleCaloClusterer &hcal)
T & operator[](int icell)
const Grid * getGrid(const std::string &type)
GridData(const Grid &grid)
Definition: CaloClusterer.h:99
bool valid_ieta_iphi(int ieta, int iphi) const
Definition: CaloClusterer.h:63
int ifind_cell(int ieta, int iphi) const
Definition: CaloClusterer.h:56
std::unique_ptr< l1t::PFClusterCollection > fetch() const
int find_cell(float eta, float phi) const override
std::vector< int > cell_map_
Definition: CaloClusterer.h:61
int neighbour(int icell, unsigned int idx) const
Definition: CaloClusterer.h:27
std::vector< float > phiWidth_
Definition: CaloClusterer.h:46
static const int phase1_nEta_
Definition: CaloClusterer.h:82
float phi(int icell) const
Definition: CaloClusterer.h:29
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static const int phase2_ietaCoarse_
Definition: CaloClusterer.h:91
std::unique_ptr< l1t::PFClusterCollection > fetchCells(bool unclusteredOnly=false, float ptMin=0.) const
static const int phase2_nEta_
Definition: CaloClusterer.h:91
std::vector< float > phi_
Definition: CaloClusterer.h:46
std::vector< std::pair< int, float > > constituents
const T & operator()(float eta, float phi) const
void add(float pt, float eta, float phi)
void correct(const Corrector &corrector)
const SingleCaloClusterer & ecal_
const T & neigh(int icell, unsigned int idx) const
static const int phase2_ietaVeryCoarse_
Definition: CaloClusterer.h:91
unsigned int size() const
void correct(const Corrector &corrector)
void grow()
possibly grow clusters by adding unclustered energy on the sides
EtGrid & raw()
non-const access to the energy: be careful to use it only before &#39;run()&#39;
float etaWidth(int icell) const
Definition: CaloClusterer.h:30
static const int phase1_ietaCoarse_
Definition: CaloClusterer.h:82
std::vector< T > data_
const EtGrid & raw() const
static const int phase2_nPhi_
Definition: CaloClusterer.h:91
float eta(int icell) const
std::unique_ptr< l1t::PFClusterCollection > fetch(float ptMin=0.) const
float phi(int icell) const
GridData< T > & operator+=(const GridData< T > &other)
SingleCaloClusterer combClusterer_
HLT enums.
GridData< int > IndexGrid
std::vector< int > ieta_
Definition: CaloClusterer.h:47
virtual int find_cell(float eta, float phi) const =0
std::vector< int > iphi_
Definition: CaloClusterer.h:47
GridData< T > & operator=(const GridData< T > &other)
Grid(unsigned int size)
Definition: CaloClusterer.h:36
const Cluster & cluster(int i) const
unsigned int size() const
Definition: CaloClusterer.h:25
std::vector< std::array< int, 8 > > neighbours_
Definition: CaloClusterer.h:48
long double T
static const float phase1_towerEtas_[phase1_nEta_]
Definition: CaloClusterer.h:83
static const float phase2_towerEtas_[phase2_nEta_]
Definition: CaloClusterer.h:92
const std::vector< Cluster > & clusters() const
float phiWidth(int icell) const
Definition: CaloClusterer.h:31
static const int phase1_nPhi_
Definition: CaloClusterer.h:82