CMS 3D CMS Logo

CaloClusterer.cc
Go to the documentation of this file.
2 
3 #include <cassert>
4 
9 
11  0, 0.087, 0.174, 0.261, 0.348, 0.435, 0.522, 0.609, 0.696, 0.783, 0.870, 0.957, 1.044, 1.131,
12  1.218, 1.305, 1.392, 1.479, 1.566, 1.653, 1.740, 1.830, 1.930, 2.043, 2.172, 2.322, 2.5, 2.650,
13  2.853, 3.139, 3.314, 3.489, 3.664, 3.839, 4.013, 4.191, 4.363, 4.538, 4.716, 4.889, 5.191};
15  0, 0.087, 0.174, 0.261, 0.348, 0.435, 0.522, 0.609, 0.696, 0.783, 0.870, 0.957, 1.044, 1.131, 1.218, 1.305,
16  1.392, 1.479, 1.564, 1.648, 1.732, 1.817, 1.901, 1.986, 2.071, 2.155, 2.240, 2.324, 2.409, 2.493, 2.577, 2.662,
17  2.747, 2.831, 2.915, 3.0, 3.139, 3.314, 3.489, 3.664, 3.839, 4.013, 4.191, 4.363, 4.538, 4.716, 4.889, 5.191};
18 
20  int nEta, int nPhi, int ietaCoarse, int ietaVeryCoarse, const float *towerEtas)
21  : Grid(2 * ((ietaCoarse - 1) * nPhi + (ietaVeryCoarse - ietaCoarse) * (nPhi / 2) +
22  (nEta - ietaVeryCoarse + 1) * (nPhi / 4))),
23  nEta_(nEta),
24  nPhi_(nPhi),
25  ietaCoarse_(ietaCoarse),
26  ietaVeryCoarse_(ietaVeryCoarse),
27  towerEtas_(towerEtas),
28  cell_map_(2 * nEta * nPhi, -1) {
29  int icell = 0;
30  for (int ie = -nEta_; ie <= nEta_; ++ie) {
31  int absie = std::abs(ie);
32  for (int iph = 1; iph <= nPhi_; ++iph) {
33  if (!valid_ieta_iphi(ie, iph))
34  continue;
35  ieta_[icell] = ie;
36  iphi_[icell] = iph;
37  float etaLo = (absie < nEta_ ? towerEtas_[absie - 1] : towerEtas_[absie - 2]);
38  float etaHi = (absie < nEta_ ? towerEtas_[absie] : towerEtas_[absie - 1]);
39  eta_[icell] = (ie > 0 ? 0.5 : -0.5) * (etaLo + etaHi);
40  etaWidth_[icell] = (etaHi - etaLo);
41  phiWidth_[icell] = 2 * M_PI / nPhi_;
42  if (absie >= ietaVeryCoarse_)
43  phiWidth_[icell] *= 4;
44  else if (absie >= ietaCoarse_)
45  phiWidth_[icell] *= 2;
46  phi_[icell] = (iph - 1) * 2 * M_PI / nPhi_ + 0.5 * phiWidth_[icell];
47  if (phi_[icell] > M_PI)
48  phi_[icell] -= 2 * M_PI;
49  std::fill(neighbours_[icell].begin(), neighbours_[icell].end(), -1);
50  cell_map_[(ie + nEta_) + 2 * nEta_ * (iph - 1)] = icell;
51  icell++;
52  }
53  }
54  assert(unsigned(icell) == ncells_);
55  // now link the cells
56  for (icell = 0; icell < int(ncells_); ++icell) {
57  int ie = ieta_[icell], iph = iphi_[icell];
58  int ineigh = 0;
59  for (int deta = -1; deta <= +1; ++deta) {
60  for (int dphi = -1; dphi <= +1; ++dphi) {
61  if (deta == 0 && dphi == 0)
62  continue;
63  neighbours_[icell][ineigh++] = imove(ie, iph, deta, dphi);
64  }
65  }
66  }
69  //for (float teta = 0; teta <= 5.0; teta += 0.02) {
70  // for (float tphi = -M_PI; tphi <= M_PI; tphi += 0.02) {
71  // find_cell(+teta, tphi);
72  // find_cell(-teta, tphi);
73  // }
74  //}
75 }
76 
77 int l1tpf_calo::Phase1GridBase::find_cell(float eta, float phi) const {
78  int ieta =
79  (eta != 0) ? std::distance(towerEtas_, std::lower_bound(towerEtas_, towerEtas_ + nEta_, std::abs(eta))) : 1;
80  if (ieta == nEta_)
81  return -1; // outside bounds
82  assert(ieta > 0 && ieta < nEta_);
83  if (ieta > nEta_)
84  ieta = nEta_;
85  if (eta < 0)
86  ieta = -ieta;
87  phi = reco::reduceRange(phi); // [-PI, PI]
88  if (phi < 0) // then bring to [0, 2*PI]
89  phi += 2 * M_PI;
90  int iphi = std::floor(phi * nPhi_ / (2 * M_PI));
91  if (phi >= 2 * M_PI)
92  iphi = nPhi_ - 1; // fix corner case due to roundings etc
93  assert(iphi < nPhi_);
94  if (std::abs(ieta) >= ietaVeryCoarse_)
95  iphi -= (iphi % 4);
96  else if (std::abs(ieta) >= ietaCoarse_)
97  iphi -= (iphi % 2);
98  iphi += 1;
100  //if (!valid_ieta_iphi(ieta,iphi)) {
101  // printf("Error in finding cell for eta %+7.4f phi %+7.4f, got ieta = %+3d iphi %2d which is not valid\n",
102  // eta, phi, ieta, iphi);
103  //}
104  assert(valid_ieta_iphi(ieta, iphi));
105  int icell = ifind_cell(ieta, iphi);
106  assert(icell != -1);
107 
109  //if (std::abs(eta - eta_[icell]) > 0.501*etaWidth_[icell] || std::abs(deltaPhi(phi, phi_[icell])) > 0.501*phiWidth_[icell]) {
110  // printf("Mismatch in finding cell for eta %+7.4f phi %+7.4f, got ieta = %+3d iphi %2d which has eta %+7.4f +- %.4f phi %+7.4f +- %.4f ; deta = %+7.4f dphi = %+7.4f\n",
111  // eta, phi, ieta, iphi, eta_[icell], etaWidth_[icell], phi_[icell], phiWidth_[icell], eta - eta_[icell], deltaPhi(phi, phi_[icell]));
112  //}
113  //assert(std::abs(eta - eta_[icell]) <= 0.5*etaWidth_[icell]);
114  //assert(std::abs(deltaPhi(phi, phi_[icell])) <= 0.5*phiWidth_[icell]);
115  return icell;
116 }
117 
118 int l1tpf_calo::Phase1GridBase::imove(int ieta, int iphi, int deta, int dphi) {
119  int ie = ieta, iph = iphi;
120  switch (deta) {
121  case -1:
122  ie = (ie == -nEta_ ? 0 : (ie == +1 ? -1 : ie - 1));
123  break;
124  case +1:
125  ie = (ie == +nEta_ ? 0 : (ie == -1 ? +1 : ie + 1));
126  break;
127  case 0:
128  break;
129  default:
130  assert(false);
131  };
132  if (ie == 0)
133  return -1;
134  switch (dphi) {
135  case -1:
136  iph = (iph == 1 ? nPhi_ : iph - 1);
137  break;
138  case +1:
139  iph = (iph == nPhi_ ? 1 : iph + 1);
140  break;
141  case 0:
142  break;
143  default:
144  assert(false);
145  };
146  if (!valid_ieta_iphi(ie, iph))
147  return -1;
148  int icell = ifind_cell(ie, iph);
149  assert(!(ie == ieta && iph == iphi));
150  assert(icell != -1);
151  assert(icell != ifind_cell(ieta, iphi));
152  return icell;
153 }
154 
156  static const Phase1Grid _phase1Grid;
157  static const Phase2Grid _phase2Grid;
158  if (type == "phase1")
159  return &_phase1Grid;
160  else if (type == "phase2")
161  return &_phase2Grid;
162  else
163  throw cms::Exception("Configuration") << "Unsupported grid type '" << type << "'\n";
164 }
165 
167  : grid_(getGrid(pset.getParameter<std::string>("grid"))),
168  rawet_(*grid_),
169  unclustered_(*grid_),
170  precluster_(*grid_),
171  clusterIndex_(*grid_),
172  cellKey_(*grid_),
173  clusters_(),
174  nullCluster_(),
175  zsEt_(pset.getParameter<double>("zsEt")),
176  seedEt_(pset.getParameter<double>("seedEt")),
177  minClusterEt_(pset.getParameter<double>("minClusterEt")),
178  minEtToGrow_(pset.existsAs<double>("minEtToGrow") ? pset.getParameter<double>("minEtToGrow") : -1),
179  energyWeightedPosition_(pset.getParameter<bool>("energyWeightedPosition")) {
180  std::string energyShareAlgo = pset.getParameter<std::string>("energyShareAlgo");
181  if (energyShareAlgo == "fractions")
183  else if (energyShareAlgo == "none")
185  else if (energyShareAlgo == "greedy")
187  else if (energyShareAlgo == "crude")
189  else
190  throw cms::Exception("Configuration") << "Unsupported energyShareAlgo '" << energyShareAlgo << "'\n";
191 }
192 
194 
196  rawet_.zero();
197  clusters_.clear();
198  clusterIndex_.fill(-1);
199 }
200 
202  unsigned int i, ncells = grid_->size();
203 
204  // kill zeros. count non-zeros, for linking later
205  cellKey_.fill(-1);
206  int key = 0;
207  for (i = 0; i < ncells; ++i) {
208  if (rawet_[i] < zsEt_) {
209  rawet_[i] = 0;
210  } else {
211  cellKey_[i] = key++;
212  }
213  }
214 
215  precluster_.clear();
216  // pre-cluster step 1: at each cell, set the value equal to itself if it's a local maxima, zero otherwise
217  // can be done in parallel on all cells
218  for (i = 0; i < ncells; ++i) {
219  if (rawet_[i] > seedEt_) {
220  precluster_[i].ptLocalMax = rawet_[i];
222  //printf(" candidate precluster pt %7.2f at %4d (ieta %+3d iphi %2d)\n", rawet_[i], i, grid_->ieta(i), grid_->iphi(i));
223  for (int ineigh = 0; ineigh <= 3; ++ineigh) {
224  if (rawet_.neigh(i, ineigh) > rawet_[i])
225  precluster_[i].ptLocalMax = 0;
227  //int ncell = grid_->neighbour(i,ineigh);
228  //if (ncell == -1) printf(" \t neigh %d is null\n", ineigh);
229  //else printf(" \t neigh %d at %4d (ieta %+3d iphi %2d) has pt %7.2f: comparison %1d \n", ineigh, ncell, grid_->ieta(ncell), grid_->iphi(ncell), rawet_[ncell], precluster_[i].ptLocalMax > 0);
230  }
231  for (int ineigh = 4; ineigh < 8; ++ineigh) {
232  if (rawet_.neigh(i, ineigh) >= rawet_[i])
233  precluster_[i].ptLocalMax = 0;
235  //int ncell = grid_->neighbour(i,ineigh);
236  //if (ncell == -1) printf(" \t neigh %d is null\n", ineigh);
237  //else printf(" \t neigh %d at %4d (ieta %+3d iphi %2d) has pt %7.2f: comparison %1d \n", ineigh, ncell, grid_->ieta(ncell), grid_->iphi(ncell), rawet_[ncell], precluster_[i].ptLocalMax > 0);
238  }
239  }
240  }
241  // pre-cluster step 2: compute information from neighbouring local max, for energy sharing purposes
242  for (i = 0; i < ncells; ++i) {
243  if (precluster_[i].ptLocalMax == 0) {
244  switch (energyShareAlgo_) {
245  case EnergyShareAlgo::Fractions: {
246  float tot = 0;
247  for (int ineigh = 0; ineigh < 8; ++ineigh) {
248  tot += precluster_.neigh(i, ineigh).ptLocalMax;
249  }
250  precluster_[i].ptOverNeighLocalMaxSum = tot ? rawet_[i] / tot : 0;
251  } break;
253  precluster_[i].ptOverNeighLocalMaxSum = rawet_[i];
254  break;
255  case EnergyShareAlgo::Greedy: {
256  float maxet = 0;
257  for (int ineigh = 0; ineigh < 8; ++ineigh) {
258  maxet = std::max(maxet, precluster_.neigh(i, ineigh).ptLocalMax);
259  }
260  precluster_[i].ptOverNeighLocalMaxSum = maxet;
261  } break;
262  case EnergyShareAlgo::Crude: {
263  int number = 0;
264  for (int ineigh = 0; ineigh < 8; ++ineigh) {
265  number += (precluster_.neigh(i, ineigh).ptLocalMax > 0);
266  }
267  precluster_[i].ptOverNeighLocalMaxSum = (number > 1 ? 0.5 : 1.0) * rawet_[i];
268  } break;
269  }
270  }
271  }
272 
273  clusterIndex_.fill(-1);
274  clusters_.clear();
275  unclustered_ = rawet_;
276  // cluster: at each localMax cell, take itself plus the weighted contributions of the neighbours
277  Cluster cluster;
278  for (i = 0; i < ncells; ++i) {
279  if (precluster_[i].ptLocalMax > 0) {
280  float myet = rawet_[i];
281  float tot = myet;
282  float avg_eta = 0;
283  float avg_phi = 0;
284  cluster.clear();
285  cluster.constituents.emplace_back(i, 1.0);
286  for (int ineigh = 0; ineigh < 8; ++ineigh) {
287  int ineighcell = grid_->neighbour(i, ineigh);
288  if (ineighcell == -1)
289  continue; // skip dummy cells
290  float fracet = 0;
291  switch (energyShareAlgo_) {
292  case EnergyShareAlgo::Fractions:
293  fracet = myet * precluster_.neigh(i, ineigh).ptOverNeighLocalMaxSum;
294  break;
296  fracet = precluster_.neigh(i, ineigh).ptOverNeighLocalMaxSum;
297  break;
298  case EnergyShareAlgo::Greedy:
299  fracet = (myet == precluster_.neigh(i, ineigh).ptOverNeighLocalMaxSum ? rawet_.neigh(i, ineigh) : 0);
300  break;
301  case EnergyShareAlgo::Crude:
302  fracet = precluster_.neigh(i, ineigh).ptOverNeighLocalMaxSum;
303  break;
304  }
305  if (fracet == 0)
306  continue;
307  tot += fracet;
308  cluster.constituents.emplace_back(ineighcell, fracet / rawet_.neigh(i, ineigh));
309  if (energyWeightedPosition_) {
310  avg_eta += fracet * (grid_->eta(ineighcell) - grid_->eta(i));
311  avg_phi += fracet * deltaPhi(grid_->phi(ineighcell), grid_->phi(i));
312  }
313  }
314  if (tot > minClusterEt_) {
315  cluster.et = tot;
316  unclustered_[i] = 0;
317  for (int ineigh = 0; ineigh < 8; ++ineigh) {
318  int ineighcell = grid_->neighbour(i, ineigh);
319  if (ineighcell == -1)
320  continue; // skip dummy cells
321  unclustered_[ineighcell] = 0;
322  }
323  if (energyWeightedPosition_) {
324  cluster.eta = grid_->eta(i) + avg_eta / tot;
325  cluster.phi = grid_->phi(i) + avg_phi / tot;
326  // wrap around phi
327  cluster.phi = reco::reduceRange(cluster.phi);
328  } else {
329  cluster.eta = grid_->eta(i);
330  cluster.phi = grid_->phi(i);
331  }
332  clusterIndex_[i] = clusters_.size();
333  clusters_.push_back(cluster);
334  }
335  }
336  }
337  if (minEtToGrow_ > 0)
338  grow();
339 }
340 
342  int selneighs[4] = {1, 3, 4, 6}; // -eta, -phi, +phi, +eta
343  std::vector<int> toreset;
344  for (Cluster &cluster : clusters_) {
345  if (cluster.et > minEtToGrow_) {
346  int i = cluster.constituents.front().first;
347  for (int side = 0; side < 4; ++side) {
348  int neigh = grid_->neighbour(i, selneighs[side]);
349  if (neigh == -1)
350  continue;
351  for (int in = 0; in < 8; ++in) {
352  int n2 = grid_->neighbour(neigh, in);
353  if (n2 == -1)
354  continue;
355  cluster.et += unclustered_[n2];
356  if (unclustered_[n2]) {
357  cluster.constituents.emplace_back(n2, 1.0);
358  toreset.push_back(n2);
359  }
360  }
361  }
362  }
363  }
364  for (int i : toreset)
365  unclustered_[i] = 0;
366 }
367 
368 std::unique_ptr<l1t::PFClusterCollection> l1tpf_calo::SingleCaloClusterer::fetchCells(bool unclusteredOnly,
369  float ptMin) const {
370  auto ret = std::make_unique<l1t::PFClusterCollection>();
371  const EtGrid &src = (unclusteredOnly ? unclustered_ : rawet_);
372  for (unsigned int i = 0, ncells = grid_->size(); i < ncells; ++i) {
373  if (src[i] <= ptMin)
374  continue;
375  if ((unclusteredOnly == false) && (ptMin == 0)) {
376  assert(cellKey_[i] == int(ret->size()));
377  }
378  ret->emplace_back(src[i], grid_->eta(i), grid_->phi(i));
379  ret->back().setHwEta(grid_->ieta(i));
380  ret->back().setHwPhi(grid_->iphi(i));
381  }
382  return ret;
383 }
384 
385 std::unique_ptr<l1t::PFClusterCollection> l1tpf_calo::SingleCaloClusterer::fetch(float ptMin) const {
386  auto ret = std::make_unique<l1t::PFClusterCollection>();
387  for (const Cluster &cluster : clusters_) {
388  if (cluster.et > ptMin) {
389  ret->emplace_back(cluster.et, cluster.eta, cluster.phi);
390  }
391  }
392  return ret;
393 }
394 
395 std::unique_ptr<l1t::PFClusterCollection> l1tpf_calo::SingleCaloClusterer::fetch(
397  auto ret = std::make_unique<l1t::PFClusterCollection>();
398  for (const Cluster &cluster : clusters_) {
399  if (cluster.et > ptMin) {
400  ret->emplace_back(cluster.et, cluster.eta, cluster.phi);
401  for (const auto &pair : cluster.constituents) {
402  edm::Ptr<l1t::PFCluster> ref(cells, cellKey_[pair.first]);
403  ret->back().addConstituent(ref, pair.second);
404  }
405  }
406  }
407  return ret;
408 }
409 
411  const SingleCaloClusterer &ecal,
412  const SingleCaloClusterer &hcal)
413  : grid_(getGrid(pset.getParameter<std::string>("grid"))),
414  ecal_(ecal),
415  hcal_(hcal),
416  clusterIndex_(*grid_),
417  clusters_(),
418  hoeCut_(pset.getParameter<double>("hoeCut")),
419  minPhotonEt_(pset.getParameter<double>("minPhotonEt")),
420  minHadronRawEt_(pset.getParameter<double>("minHadronRawEt")),
421  minHadronEt_(pset.getParameter<double>("minHadronEt")),
422  noEmInHGC_(pset.getParameter<bool>("noEmInHGC")) {
423  if (grid_ != &ecal.raw().grid())
424  throw cms::Exception("LogicError", "Inconsistent grid between ecal and linker\n");
425  if (grid_ != &hcal.raw().grid())
426  throw cms::Exception("LogicError", "Inconsistent grid between hcal and linker\n");
427 }
428 
430 
431 std::unique_ptr<l1t::PFClusterCollection> l1tpf_calo::SimpleCaloLinkerBase::fetch() const {
433  return fetch(ecal, hcal);
434 }
435 
436 std::unique_ptr<l1t::PFClusterCollection> l1tpf_calo::SimpleCaloLinkerBase::fetch(
439  bool setRefs = (ecal.isValid() && hcal.isValid());
440  auto ret = std::make_unique<l1t::PFClusterCollection>();
441  for (const CombinedCluster &cluster : clusters_) {
442  if (cluster.et > 0) {
443  bool photon = (cluster.hcal_et < hoeCut_ * cluster.ecal_et);
444  if (photon && noEmInHGC_) {
445  if (std::abs(cluster.eta) > 1.5 && std::abs(cluster.eta) < 3.0) { // 1.5-3 = eta range of HGCal
446  continue;
447  }
448  }
449  if (cluster.et > (photon ? minPhotonEt_ : minHadronEt_)) {
450  ret->emplace_back(cluster.et,
451  cluster.eta,
452  cluster.phi,
453  cluster.ecal_et > 0 ? std::max(cluster.et - cluster.ecal_et, 0.f) / cluster.ecal_et : -1,
454  photon);
455  if (setRefs) {
456  for (const auto &pair : cluster.constituents) {
457  assert(pair.first != 0);
458  if (pair.first > 0) { // 1+hcal index
459  ret->back().addConstituent(edm::Ptr<l1t::PFCluster>(hcal, +pair.first - 1), pair.second);
460  } else { // -1-ecal index
461  ret->back().addConstituent(edm::Ptr<l1t::PFCluster>(ecal, -pair.first + 1), pair.second);
462  }
463  }
464  }
465  }
466  }
467  }
468  return ret;
469 }
470 
472  const SingleCaloClusterer &ecal,
473  const SingleCaloClusterer &hcal)
474  : SimpleCaloLinkerBase(pset, ecal, hcal), ecalToHCal_(*grid_) {}
475 
477 
479  clearBase();
480  ecalToHCal_.clear();
481 }
482 
484  unsigned int i, ncells = grid_->size();
485 
486  const EtGrid &hraw = hcal_.raw();
487  const IndexGrid &ecals = ecal_.indexGrid();
488  const IndexGrid &hcals = hcal_.indexGrid();
489 
490  // for each ECal cluster, get the corresponding HCal cluster and the sum of the neighbour HCal clusters
491  ecalToHCal_.clear();
492  for (i = 0; i < ncells; ++i) {
493  if (ecals[i] >= 0) {
494  if (hcals[i] >= 0) {
495  ecalToHCal_[i].ptLocalMax = hcal_.cluster(i).et;
496  } else {
497  float tot = 0;
498  for (int ineigh = 0; ineigh < 8; ++ineigh) {
499  tot += hcal_.cluster(grid_->neighbour(i, ineigh)).et;
500  }
501  ecalToHCal_[i].ptOverNeighLocalMaxSum = tot ? ecal_.cluster(i).et / tot : 0;
502  }
503  }
504  }
505 
506  clusterIndex_.fill(-1);
507  clusters_.clear();
508  CombinedCluster cluster;
509  // promote HCal clusters to final clusters
510  for (i = 0; i < ncells; ++i) {
511  if (hcals[i] >= 0) {
512  const Cluster &hcal = hcal_.cluster(i);
513  cluster.clear();
514  cluster.constituents.emplace_back(+i + 1, 1);
515  if (ecalToHCal_[i].ptLocalMax > 0) {
516  // direct linking is easy
517  const Cluster &ecal = ecal_.cluster(i);
518  if (ecal.et + hcal.et > minHadronRawEt_) {
519  cluster.ecal_et = ecal.et;
520  cluster.hcal_et = hcal.et;
521  cluster.et = cluster.ecal_et + cluster.hcal_et;
522  float wecal = cluster.ecal_et / cluster.et, whcal = 1.0 - wecal;
523  cluster.eta = ecal.eta * wecal + hcal.eta * whcal;
524  cluster.phi = ecal.phi * wecal + hcal.phi * whcal;
525  // wrap around phi
526  cluster.phi = reco::reduceRange(cluster.phi);
527  cluster.constituents.emplace_back(-i - 1, 1);
528  }
529  } else {
530  // sidewas linking is more annonying
531  float myet = hcal.et;
532  float etot = 0;
533  float avg_eta = 0;
534  float avg_phi = 0;
535  for (int ineigh = 0; ineigh < 8; ++ineigh) {
536  int ineighcell = grid_->neighbour(i, ineigh);
537  if (ineighcell == -1)
538  continue; // skip dummy cells
539  float fracet = myet * ecalToHCal_.neigh(i, ineigh).ptOverNeighLocalMaxSum;
540  if (fracet == 0)
541  continue;
542  etot += fracet;
543  avg_eta += fracet * (grid_->eta(ineighcell) - grid_->eta(i));
544  avg_phi += fracet * deltaPhi(grid_->phi(ineighcell), grid_->phi(i));
545  cluster.constituents.emplace_back(-i - 1, fracet / ecal_.cluster(ineighcell).et);
546  }
547  if (myet + etot > minHadronRawEt_) {
548  cluster.hcal_et = hcal.et;
549  cluster.ecal_et = etot;
550  cluster.et = myet + etot;
551  cluster.eta = hcal.eta + avg_eta / cluster.et;
552  cluster.phi = hcal.phi + avg_phi / cluster.et;
553  // wrap around phi
554  cluster.phi = reco::reduceRange(cluster.phi);
555  }
556  }
557  if (cluster.et > 0) {
558  clusterIndex_[i] = clusters_.size();
559  clusters_.push_back(cluster);
560  }
561  }
562  }
563 
564  // promote Unlinked ECal clusters to final clusters
565  for (i = 0; i < ncells; ++i) {
566  if (ecals[i] >= 0 && ecalToHCal_[i].ptLocalMax == 0 && ecalToHCal_[i].ptOverNeighLocalMaxSum == 0) {
567  cluster.clear();
568  const Cluster &ecal = ecal_.cluster(i);
569  cluster.ecal_et = ecal.et;
570  cluster.hcal_et = hraw[i];
571  cluster.et = cluster.ecal_et + cluster.hcal_et;
572  cluster.eta = ecal.eta;
573  cluster.phi = ecal.phi;
574  cluster.constituents.emplace_back(-i - 1, 1);
575  clusterIndex_[i] = clusters_.size();
576  clusters_.push_back(cluster);
577  }
578  }
579 }
580 
582  const SingleCaloClusterer &ecal,
583  const SingleCaloClusterer &hcal)
584  : SimpleCaloLinkerBase(pset, ecal, hcal), combClusterer_(pset) {}
585 
587 
589  clearBase();
590  combClusterer_.clear();
591 }
592 
594  combClusterer_.clear();
595 
596  const EtGrid &hraw = hcal_.raw();
597  const EtGrid &eraw = ecal_.raw();
598  combClusterer_.raw() = eraw;
599  combClusterer_.raw() += hraw;
600 
601  combClusterer_.run();
602  clusterIndex_ = combClusterer_.indexGrid();
603  const std::vector<Cluster> &clustersSrc = combClusterer_.clusters();
604  unsigned int nclust = clustersSrc.size();
605  clusters_.resize(nclust);
606  for (unsigned int ic = 0; ic < nclust; ++ic) {
607  const Cluster &src = clustersSrc[ic];
608  CombinedCluster &dst = clusters_[ic];
609  dst.et = src.et;
610  dst.eta = src.eta;
611  dst.phi = src.phi;
612  dst.ecal_et = 0;
613  dst.hcal_et = 0;
614  for (const auto &pair : src.constituents) {
615  if (eraw[pair.first]) {
616  dst.ecal_et += pair.second * eraw[pair.first];
617  dst.constituents.emplace_back(-pair.first - 1, pair.second);
618  }
619  if (hraw[pair.first]) {
620  dst.hcal_et += pair.second * hraw[pair.first];
621  dst.constituents.emplace_back(+pair.first + 1, pair.second);
622  }
623  }
624  }
625 }
626 
627 std::unique_ptr<l1tpf_calo::SimpleCaloLinkerBase> l1tpf_calo::makeCaloLinker(const edm::ParameterSet &pset,
628  const SingleCaloClusterer &ecal,
629  const SingleCaloClusterer &hcal) {
630  const std::string &algo = pset.getParameter<std::string>("algo");
631  if (algo == "simple") {
632  return std::make_unique<l1tpf_calo::SimpleCaloLinker>(pset, ecal, hcal);
633  } else if (algo == "flat") {
634  return std::make_unique<l1tpf_calo::FlatCaloLinker>(pset, ecal, hcal);
635  } else {
636  throw cms::Exception("Configuration") << "Unsupported linker algo '" << algo << "'\n";
637  }
638 }
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition: runTheMatrix.py:355
l1tpf_calo::CombinedCluster::clear
void clear()
Definition: CaloClusterer.h:172
muons2muons_cfi.photon
photon
Definition: muons2muons_cfi.py:28
l1tpf_calo::Grid::neighbours_
std::vector< std::array< int, 8 > > neighbours_
Definition: CaloClusterer.h:48
l1tpf_calo::SimpleCaloLinkerBase::~SimpleCaloLinkerBase
virtual ~SimpleCaloLinkerBase()
Definition: CaloClusterer.cc:429
electrons_cff.bool
bool
Definition: electrons_cff.py:372
l1tpf_calo::GridData::clear
void clear()
Definition: CaloClusterer.h:139
mps_fire.i
i
Definition: mps_fire.py:355
l1tpf_calo::Phase1GridBase::nEta_
const int nEta_
Definition: CaloClusterer.h:59
l1tpf_calo::Phase1GridBase::towerEtas_
const float * towerEtas_
Definition: CaloClusterer.h:60
l1tpf_calo::Phase1GridBase::nPhi_
const int nPhi_
Definition: CaloClusterer.h:59
l1tpf_calo::Phase1GridBase::imove
int imove(int ieta, int iphi, int deta, int dphi)
Definition: CaloClusterer.cc:118
l1tpf_calo::getGrid
const Grid * getGrid(const std::string &type)
Definition: CaloClusterer.cc:155
hcal
Definition: ConfigurationDatabase.cc:13
l1tpf_calo::CombinedCluster::ecal_et
float ecal_et
Definition: CaloClusterer.h:171
l1tpf_calo::Grid
Definition: CaloClusterer.h:22
l1tpf_calo::Grid::ncells_
unsigned int ncells_
Definition: CaloClusterer.h:45
l1tpf_calo::FlatCaloLinker::clear
void clear() override
Definition: CaloClusterer.cc:588
ptMin
constexpr float ptMin
Definition: PhotonIDValueMapProducer.cc:153
cms::cuda::assert
assert(be >=bs)
HLT_2018_cff.distance
distance
Definition: HLT_2018_cff.py:6417
l1tpf_calo::Grid::phi_
std::vector< float > phi_
Definition: CaloClusterer.h:46
l1tpf_calo::Phase2Grid::phase2_nEta_
static const int phase2_nEta_
Definition: CaloClusterer.h:91
l1tpf_calo::Phase1Grid::phase1_towerEtas_
static const float phase1_towerEtas_[phase1_nEta_]
Definition: CaloClusterer.h:83
l1tpf_calo::FlatCaloLinker::run
void run() override
Definition: CaloClusterer.cc:593
pfClustersFromCombinedCalo_cfi.energyShareAlgo
energyShareAlgo
Definition: pfClustersFromCombinedCalo_cfi.py:22
l1tpf_calo::SimpleCaloLinker::SimpleCaloLinker
SimpleCaloLinker(const edm::ParameterSet &pset, const SingleCaloClusterer &ecal, const SingleCaloClusterer &hcal)
Definition: CaloClusterer.cc:471
l1tpf_calo::makeCaloLinker
std::unique_ptr< SimpleCaloLinkerBase > makeCaloLinker(const edm::ParameterSet &pset, const SingleCaloClusterer &ecal, const SingleCaloClusterer &hcal)
Definition: CaloClusterer.cc:627
end
#define end
Definition: vmac.h:39
l1tpf_calo::Cluster::constituents
std::vector< std::pair< int, float > > constituents
Definition: CaloClusterer.h:163
cond::persistency::fetch
std::pair< std::string, std::shared_ptr< void > > fetch(const cond::Hash &payloadId, Session &session)
Definition: CondDBFetch.cc:336
l1tpf_calo::SimpleCaloLinkerBase::SimpleCaloLinkerBase
SimpleCaloLinkerBase(const edm::ParameterSet &pset, const SingleCaloClusterer &ecal, const SingleCaloClusterer &hcal)
Definition: CaloClusterer.cc:410
LEDCalibrationChannels.iphi
iphi
Definition: LEDCalibrationChannels.py:64
patCandidatesForDimuonsSequences_cff.hcal
hcal
Definition: patCandidatesForDimuonsSequences_cff.py:37
cmsdt::algo
algo
Definition: constants.h:164
l1tpf_calo::Cluster
Definition: CaloClusterer.h:160
None
Definition: APVGainStruct.h:52
l1tpf_calo::Cluster::clear
void clear()
Definition: CaloClusterer.h:164
contentValuesFiles.number
number
Definition: contentValuesFiles.py:53
l1tpf_calo::Grid::etaWidth_
std::vector< float > etaWidth_
Definition: CaloClusterer.h:46
l1tpf_calo::SingleCaloClusterer::grow
void grow()
possibly grow clusters by adding unclustered energy on the sides
Definition: CaloClusterer.cc:341
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
PVValHelper::eta
Definition: PVValidationHelpers.h:69
l1tpf_calo::FlatCaloLinker::FlatCaloLinker
FlatCaloLinker(const edm::ParameterSet &pset, const SingleCaloClusterer &ecal, const SingleCaloClusterer &hcal)
Definition: CaloClusterer.cc:581
l1tpf_calo::Phase1GridBase::cell_map_
std::vector< int > cell_map_
Definition: CaloClusterer.h:61
l1tpf_calo::SimpleCaloLinkerBase::fetch
std::unique_ptr< l1t::PFClusterCollection > fetch() const
Definition: CaloClusterer.cc:431
l1tpf_calo::Phase1GridBase::valid_ieta_iphi
bool valid_ieta_iphi(int ieta, int iphi) const
Definition: CaloClusterer.h:63
l1tpf_calo::SingleCaloClusterer::fetchCells
std::unique_ptr< l1t::PFClusterCollection > fetchCells(bool unclusteredOnly=false, float ptMin=0.) const
Definition: CaloClusterer.cc:368
l1tpf_calo::Grid::ieta_
std::vector< int > ieta_
Definition: CaloClusterer.h:47
cuda_std::lower_bound
__host__ constexpr __device__ RandomIt lower_bound(RandomIt first, RandomIt last, const T &value, Compare comp={})
Definition: cudastdAlgorithm.h:27
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
l1tpf_calo::Phase1GridBase::find_cell
int find_cell(float eta, float phi) const override
Definition: CaloClusterer.cc:77
l1tpf_calo::SingleCaloClusterer::~SingleCaloClusterer
~SingleCaloClusterer()
Definition: CaloClusterer.cc:193
ntuplemaker.fill
fill
Definition: ntuplemaker.py:304
l1tpf_calo::Phase1GridBase::ietaVeryCoarse_
const int ietaVeryCoarse_
Definition: CaloClusterer.h:59
LEDCalibrationChannels.ieta
ieta
Definition: LEDCalibrationChannels.py:63
edm::ParameterSet
Definition: ParameterSet.h:36
TrackRefitter_38T_cff.src
src
Definition: TrackRefitter_38T_cff.py:24
l1tpf_calo::SingleCaloClusterer::clear
void clear()
Definition: CaloClusterer.cc:195
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
l1tpf_calo::CombinedCluster
Definition: CaloClusterer.h:170
recoMuon::in
Definition: RecoMuonEnumerators.h:6
createfilelist.int
int
Definition: createfilelist.py:10
RefToPtr.h
l1tpf_calo::Grid::phiWidth_
std::vector< float > phiWidth_
Definition: CaloClusterer.h:46
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:50
HLT_2018_cff.nEta
nEta
Definition: HLT_2018_cff.py:5271
l1tpf_calo::Cluster::eta
float eta
Definition: CaloClusterer.h:162
EgHLTOffHistBins_cfi.et
et
Definition: EgHLTOffHistBins_cfi.py:8
l1tpf_calo::Phase2Grid::phase2_towerEtas_
static const float phase2_towerEtas_[phase2_nEta_]
Definition: CaloClusterer.h:92
l1tpf_calo::SingleCaloClusterer::EnergyShareAlgo::Fractions
l1tpf_calo::Phase1GridBase::ietaCoarse_
const int ietaCoarse_
Definition: CaloClusterer.h:59
l1tpf_calo::SingleCaloClusterer::run
void run()
Definition: CaloClusterer.cc:201
edm::Ptr
Definition: AssociationVector.h:31
l1tpf_calo::Cluster::et
float et
Definition: CaloClusterer.h:162
l1tpf_calo::SimpleCaloLinkerBase::grid_
const Grid * grid_
Definition: CaloClusterer.h:261
l1tpf_calo::GridData< float >
l1tpf_calo::Phase2Grid
Definition: CaloClusterer.h:85
l1tpf_calo::SimpleCaloLinker::run
void run() override
Definition: CaloClusterer.cc:483
l1tpf_calo::Grid::eta_
std::vector< float > eta_
Definition: CaloClusterer.h:46
type
type
Definition: HCALResponse.h:21
l1tpf_calo::SingleCaloClusterer::EnergyShareAlgo::Greedy
std
Definition: JetResolutionObject.h:76
edm::OrphanHandle
Definition: EDProductfwd.h:39
l1tpf_calo::Cluster::phi
float phi
Definition: CaloClusterer.h:162
l1tpf_calo::SingleCaloClusterer::EnergyShareAlgo::Crude
Exception
Definition: hltDiff.cc:246
postprocess-scan-build.cells
cells
Definition: postprocess-scan-build.py:13
l1tpf_calo::SingleCaloClusterer::EnergyShareAlgo::None
l1tpf_calo::SimpleCaloLinkerBase
Definition: CaloClusterer.h:235
l1tpf_calo::FlatCaloLinker::~FlatCaloLinker
~FlatCaloLinker() override
Definition: CaloClusterer.cc:586
DeadROCCounter.nclust
nclust
Definition: DeadROCCounter.py:66
Exception.h
l1tpf_calo::SingleCaloClusterer
Definition: CaloClusterer.h:180
bsc_activity_cfg.ecal
ecal
Definition: bsc_activity_cfg.py:25
l1tpf_calo::Grid::iphi_
std::vector< int > iphi_
Definition: CaloClusterer.h:47
HLT_2018_cff.nPhi
nPhi
Definition: HLT_2018_cff.py:5272
cms::Exception
Definition: Exception.h:70
math::cholesky::dst
M2 & dst
Definition: choleskyInversion.h:158
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
crabWrapper.key
key
Definition: crabWrapper.py:19
l1tpf_calo::CombinedCluster::hcal_et
float hcal_et
Definition: CaloClusterer.h:171
l1tpf_calo::SimpleCaloLinker::~SimpleCaloLinker
~SimpleCaloLinker() override
Definition: CaloClusterer.cc:476
l1tpf_calo::SingleCaloClusterer::energyShareAlgo_
EnergyShareAlgo energyShareAlgo_
Definition: CaloClusterer.h:231
l1tpf_calo::SingleCaloClusterer::fetch
std::unique_ptr< l1t::PFClusterCollection > fetch(float ptMin=0.) const
Definition: CaloClusterer.cc:385
l1tpf_calo::SingleCaloClusterer::SingleCaloClusterer
SingleCaloClusterer(const edm::ParameterSet &pset)
Definition: CaloClusterer.cc:166
begin
#define begin
Definition: vmac.h:32
l1tpf_calo::Phase1GridBase::Phase1GridBase
Phase1GridBase(int nEta, int nPhi, int ietaCoarse, int ietaVeryCoarse, const float *towerEtas)
Definition: CaloClusterer.cc:19
l1tpf_calo::SimpleCaloLinker::clear
void clear() override
Definition: CaloClusterer.cc:478
l1tpf_calo::Phase1Grid
Definition: CaloClusterer.h:76
HiEvtPlane_cfi.maxet
maxet
Definition: HiEvtPlane_cfi.py:18
deltaPhi.h
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
CaloClusterer.h
reco::reduceRange
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
l1tpf_calo::Phase1Grid::phase1_nEta_
static const int phase1_nEta_
Definition: CaloClusterer.h:82