CMS 3D CMS Logo

HGCalCLUEAlgo.cc
Go to the documentation of this file.
2 
3 // Geometry
8 
10 //
12 #include "tbb/task_arena.h"
13 #include "tbb/tbb.h"
14 #include <limits>
15 
16 using namespace hgcal_clustering;
17 
18 template <typename T>
20  cells_.clear();
21  numberOfClustersPerLayer_.clear();
22  cells_.resize(2 * (maxlayer_ + 1));
23  numberOfClustersPerLayer_.resize(2 * (maxlayer_ + 1), 0);
24 }
25 
26 template <typename T>
28  // loop over all hits and create the Hexel structure, skip energies below ecut
29 
30  if (dependSensor_) {
31  // for each layer and wafer calculate the thresholds (sigmaNoise and energy)
32  // once
33  computeThreshold();
34  }
35 
36  for (unsigned int i = 0; i < hits.size(); ++i) {
37  const HGCRecHit& hgrh = hits[i];
38  DetId detid = hgrh.detid();
39  unsigned int layerOnSide = (rhtools_.getLayerWithOffset(detid) - 1);
40 
41  // set sigmaNoise default value 1 to use kappa value directly in case of
42  // sensor-independent thresholds
43  float sigmaNoise = 1.f;
44  if (dependSensor_) {
45  int thickness_index = rhtools_.getSiThickIndex(detid);
46  if (thickness_index == -1)
47  thickness_index = maxNumberOfThickIndices_;
48 
49  double storedThreshold = thresholds_[layerOnSide][thickness_index];
50  if (detid.det() == DetId::HGCalHSi || detid.subdetId() == HGCHEF) {
51  storedThreshold = thresholds_[layerOnSide][thickness_index + deltasi_index_regemfac_];
52  }
53  sigmaNoise = v_sigmaNoise_[layerOnSide][thickness_index];
54 
55  if (hgrh.energy() < storedThreshold)
56  continue; // this sets the ZS threshold at ecut times the sigma noise
57  // for the sensor
58  }
59  if (!dependSensor_ && hgrh.energy() < ecut_)
60  continue;
61  const GlobalPoint position(rhtools_.getPosition(detid));
62  int offset = ((rhtools_.zside(detid) + 1) >> 1) * maxlayer_;
63  int layer = layerOnSide + offset;
64 
65  cells_[layer].detid.emplace_back(detid);
66  cells_[layer].x.emplace_back(position.x());
67  cells_[layer].y.emplace_back(position.y());
68  if (!rhtools_.isOnlySilicon(layer)) {
69  cells_[layer].isSi.emplace_back(rhtools_.isSilicon(detid));
70  cells_[layer].eta.emplace_back(position.eta());
71  cells_[layer].phi.emplace_back(position.phi());
72  } // else, isSilicon == true and eta phi values will not be used
73  cells_[layer].weight.emplace_back(hgrh.energy());
74  cells_[layer].sigmaNoise.emplace_back(sigmaNoise);
75  }
76 }
77 
78 template <typename T>
80  auto cellsSize = cells_[l].detid.size();
81  cells_[l].rho.resize(cellsSize, 0.f);
82  cells_[l].delta.resize(cellsSize, 9999999);
83  cells_[l].nearestHigher.resize(cellsSize, -1);
84  cells_[l].clusterIndex.resize(cellsSize, -1);
85  cells_[l].followers.resize(cellsSize);
86  cells_[l].isSeed.resize(cellsSize, false);
87  if (rhtools_.isOnlySilicon(l)) {
88  cells_[l].isSi.resize(cellsSize, true);
89  cells_[l].eta.resize(cellsSize, 0.f);
90  cells_[l].phi.resize(cellsSize, 0.f);
91  }
92 }
93 
94 // Create a vector of Hexels associated to one cluster from a collection of
95 // HGCalRecHits - this can be used directly to make the final cluster list -
96 // this method can be invoked multiple times for the same event with different
97 // input (reset should be called between events)
98 template <typename T>
100  // assign all hits in each layer to a cluster core
101  tbb::this_task_arena::isolate([&] {
102  tbb::parallel_for(size_t(0), size_t(2 * maxlayer_ + 2), [&](size_t i) {
103  prepareDataStructures(i);
104  T lt;
105  lt.clear();
106  lt.fill(cells_[i].x, cells_[i].y, cells_[i].eta, cells_[i].phi, cells_[i].isSi);
107  float delta_c; // maximum search distance (critical distance) for local
108  // density calculation
109  if (i % maxlayer_ < lastLayerEE_)
110  delta_c = vecDeltas_[0];
111  else if (i % maxlayer_ < (firstLayerBH_ - 1))
112  delta_c = vecDeltas_[1];
113  else
114  delta_c = vecDeltas_[2];
115  float delta_r = vecDeltas_[3];
116  LogDebug("HGCalCLUEAlgo") << "maxlayer: " << maxlayer_ << " lastLayerEE: " << lastLayerEE_
117  << " firstLayerBH: " << firstLayerBH_ << "\n";
118 
119  calculateLocalDensity(lt, i, delta_c, delta_r);
120  calculateDistanceToHigher(lt, i, delta_c, delta_r);
121  numberOfClustersPerLayer_[i] = findAndAssignClusters(i, delta_c, delta_r);
122  });
123  });
124  //Now that we have the density per point we can store it
125  for (unsigned int i = 0; i < 2 * maxlayer_ + 2; ++i) {
126  setDensity(i);
127  }
128 }
129 
130 template <typename T>
131 std::vector<reco::BasicCluster> HGCalCLUEAlgoT<T>::getClusters(bool) {
132  std::vector<int> offsets(numberOfClustersPerLayer_.size(), 0);
133 
134  int maxClustersOnLayer = numberOfClustersPerLayer_[0];
135 
136  for (unsigned layerId = 1; layerId < offsets.size(); ++layerId) {
137  offsets[layerId] = offsets[layerId - 1] + numberOfClustersPerLayer_[layerId - 1];
138 
139  maxClustersOnLayer = std::max(maxClustersOnLayer, numberOfClustersPerLayer_[layerId]);
140  }
141 
142  auto totalNumberOfClusters = offsets.back() + numberOfClustersPerLayer_.back();
143  clusters_v_.resize(totalNumberOfClusters);
144  std::vector<std::vector<int>> cellsIdInCluster;
145  cellsIdInCluster.reserve(maxClustersOnLayer);
146 
147  for (unsigned int layerId = 0; layerId < 2 * maxlayer_ + 2; ++layerId) {
148  cellsIdInCluster.resize(numberOfClustersPerLayer_[layerId]);
149  auto& cellsOnLayer = cells_[layerId];
150  unsigned int numberOfCells = cellsOnLayer.detid.size();
151  auto firstClusterIdx = offsets[layerId];
152 
153  for (unsigned int i = 0; i < numberOfCells; ++i) {
154  auto clusterIndex = cellsOnLayer.clusterIndex[i];
155  if (clusterIndex != -1)
156  cellsIdInCluster[clusterIndex].push_back(i);
157  }
158 
159  std::vector<std::pair<DetId, float>> thisCluster;
160 
161  for (auto& cl : cellsIdInCluster) {
162  auto position = calculatePosition(cl, layerId);
163  float energy = 0.f;
164  int seedDetId = -1;
165 
166  for (auto cellIdx : cl) {
167  energy += cellsOnLayer.weight[cellIdx];
168  thisCluster.emplace_back(cellsOnLayer.detid[cellIdx], 1.f);
169  if (cellsOnLayer.isSeed[cellIdx]) {
170  seedDetId = cellsOnLayer.detid[cellIdx];
171  }
172  }
173  auto globalClusterIndex = cellsOnLayer.clusterIndex[cl[0]] + firstClusterIdx;
174 
175  clusters_v_[globalClusterIndex] =
177  clusters_v_[globalClusterIndex].setSeed(seedDetId);
178  thisCluster.clear();
179  }
180 
181  cellsIdInCluster.clear();
182  }
183  return clusters_v_;
184 }
185 
186 template <typename T>
187 math::XYZPoint HGCalCLUEAlgoT<T>::calculatePosition(const std::vector<int>& v, const unsigned int layerId) const {
188  float total_weight = 0.f;
189  float x = 0.f;
190  float y = 0.f;
191 
192  unsigned int maxEnergyIndex = 0;
193  float maxEnergyValue = 0.f;
194 
195  auto& cellsOnLayer = cells_[layerId];
196 
197  // loop over hits in cluster candidate
198  // determining the maximum energy hit
199  for (auto i : v) {
200  total_weight += cellsOnLayer.weight[i];
201  if (cellsOnLayer.weight[i] > maxEnergyValue) {
202  maxEnergyValue = cellsOnLayer.weight[i];
203  maxEnergyIndex = i;
204  }
205  }
206 
207  // Si cell or Scintillator. Used to set approach and parameters
208  auto thick = rhtools_.getSiThickIndex(cellsOnLayer.detid[maxEnergyIndex]);
209  bool isSiliconCell = (thick != -1);
210 
211  // TODO: this is recomputing everything twice and overwriting the position with log weighting position
212  if (isSiliconCell) {
213  float total_weight_log = 0.f;
214  float x_log = 0.f;
215  float y_log = 0.f;
216  for (auto i : v) {
217  //for silicon only just use 1+6 cells = 1.3cm for all thicknesses
218  if (distance2(i, maxEnergyIndex, layerId, false) > positionDeltaRho2_)
219  continue;
220  float rhEnergy = cellsOnLayer.weight[i];
221  float Wi = std::max(thresholdW0_[thick] + std::log(rhEnergy / total_weight), 0.);
222  x_log += cellsOnLayer.x[i] * Wi;
223  y_log += cellsOnLayer.y[i] * Wi;
224  total_weight_log += Wi;
225  }
226 
227  total_weight = total_weight_log;
228  x = x_log;
229  y = y_log;
230  } else {
231  for (auto i : v) {
232  float rhEnergy = cellsOnLayer.weight[i];
233 
234  x += cellsOnLayer.x[i] * rhEnergy;
235  y += cellsOnLayer.y[i] * rhEnergy;
236  }
237  }
238  if (total_weight != 0.) {
239  float inv_tot_weight = 1.f / total_weight;
240  return math::XYZPoint(
241  x * inv_tot_weight, y * inv_tot_weight, rhtools_.getPosition(cellsOnLayer.detid[maxEnergyIndex]).z());
242  } else
243  return math::XYZPoint(0.f, 0.f, 0.f);
244 }
245 
246 template <typename T>
247 void HGCalCLUEAlgoT<T>::calculateLocalDensity(const T& lt, const unsigned int layerId, float delta_c, float delta_r) {
248  auto& cellsOnLayer = cells_[layerId];
249  unsigned int numberOfCells = cellsOnLayer.detid.size();
250  bool isOnlySi(false);
251  if (rhtools_.isOnlySilicon(layerId))
252  isOnlySi = true;
253 
254  for (unsigned int i = 0; i < numberOfCells; i++) {
255  bool isSi = isOnlySi || cellsOnLayer.isSi[i];
256  if (isSi) {
257  float delta = delta_c;
258  std::array<int, 4> search_box = lt.searchBox(
259  cellsOnLayer.x[i] - delta, cellsOnLayer.x[i] + delta, cellsOnLayer.y[i] - delta, cellsOnLayer.y[i] + delta);
260 
261  for (int xBin = search_box[0]; xBin < search_box[1] + 1; ++xBin) {
262  for (int yBin = search_box[2]; yBin < search_box[3] + 1; ++yBin) {
263  int binId = lt.getGlobalBinByBin(xBin, yBin);
264  size_t binSize = lt[binId].size();
265 
266  for (unsigned int j = 0; j < binSize; j++) {
267  unsigned int otherId = lt[binId][j];
268  bool otherSi = isOnlySi || cellsOnLayer.isSi[otherId];
269  if (otherSi) { //silicon cells cannot talk to scintillator cells
270  if (distance(i, otherId, layerId, false) < delta) {
271  cellsOnLayer.rho[i] += (i == otherId ? 1.f : 0.5f) * cellsOnLayer.weight[otherId];
272  }
273  }
274  }
275  }
276  }
277  } else {
278  float delta = delta_r;
279  std::array<int, 4> search_box = lt.searchBoxEtaPhi(cellsOnLayer.eta[i] - delta,
280  cellsOnLayer.eta[i] + delta,
281  cellsOnLayer.phi[i] - delta,
282  cellsOnLayer.phi[i] + delta);
283  cellsOnLayer.rho[i] += cellsOnLayer.weight[i];
284  float northeast(0), northwest(0), southeast(0), southwest(0), all(0);
285 
286  for (int etaBin = search_box[0]; etaBin < search_box[1] + 1; ++etaBin) {
287  for (int phiBin = search_box[2]; phiBin < search_box[3] + 1; ++phiBin) {
288  int binId = lt.getGlobalBinByBinEtaPhi(etaBin, phiBin);
289  size_t binSize = lt[binId].size();
290 
291  for (unsigned int j = 0; j < binSize; j++) {
292  unsigned int otherId = lt[binId][j];
293  if (!cellsOnLayer.isSi[otherId]) { //scintillator cells cannot talk to silicon cells
294  if (distance(i, otherId, layerId, true) < delta) {
295  int iPhi = HGCScintillatorDetId(cellsOnLayer.detid[i]).iphi();
296  int otherIPhi = HGCScintillatorDetId(cellsOnLayer.detid[otherId]).iphi();
297  int iEta = HGCScintillatorDetId(cellsOnLayer.detid[i]).ieta();
298  int otherIEta = HGCScintillatorDetId(cellsOnLayer.detid[otherId]).ieta();
299  int dIPhi = otherIPhi - iPhi;
300  dIPhi += abs(dIPhi) < 2 ? 0
301  : dIPhi < 0 ? scintMaxIphi_
302  : -scintMaxIphi_; // cells with iPhi=288 and iPhi=1 should be neiboring cells
303  int dIEta = otherIEta - iEta;
304  LogDebug("HGCalCLUEAlgo") << " Debugging calculateLocalDensity for Scintillator: \n"
305  << " cell: " << otherId << " energy: " << cellsOnLayer.weight[otherId]
306  << " otherIPhi: " << otherIPhi << " iPhi: " << iPhi
307  << " otherIEta: " << otherIEta << " iEta: " << iEta << "\n";
308 
309  if (otherId != i) {
310  auto neighborCellContribution = 0.5f * cellsOnLayer.weight[otherId];
311  all += neighborCellContribution;
312  if (dIPhi >= 0 && dIEta >= 0)
313  northeast += neighborCellContribution;
314  if (dIPhi <= 0 && dIEta >= 0)
315  southeast += neighborCellContribution;
316  if (dIPhi >= 0 && dIEta <= 0)
317  northwest += neighborCellContribution;
318  if (dIPhi <= 0 && dIEta <= 0)
319  southwest += neighborCellContribution;
320  }
321  LogDebug("HGCalCLUEAlgo") << " Debugging calculateLocalDensity for Scintillator: \n"
322  << " northeast: " << northeast << " southeast: " << southeast
323  << " northwest: " << northwest << " southwest: " << southwest << "\n";
324  }
325  }
326  }
327  }
328  }
329  float neighborsval = (std::max(northeast, northwest) > std::max(southeast, southwest))
330  ? std::max(northeast, northwest)
331  : std::max(southeast, southwest);
332  if (use2x2_)
333  cellsOnLayer.rho[i] += neighborsval;
334  else
335  cellsOnLayer.rho[i] += all;
336  }
337  LogDebug("HGCalCLUEAlgo") << "Debugging calculateLocalDensity: \n"
338  << " cell: " << i << " isSilicon: " << cellsOnLayer.isSi[i]
339  << " eta: " << cellsOnLayer.eta[i] << " phi: " << cellsOnLayer.phi[i]
340  << " energy: " << cellsOnLayer.weight[i] << " density: " << cellsOnLayer.rho[i] << "\n";
341  }
342 }
343 
344 template <typename T>
346  const unsigned int layerId,
347  float delta_c,
348  float delta_r) {
349  auto& cellsOnLayer = cells_[layerId];
350  unsigned int numberOfCells = cellsOnLayer.detid.size();
351  bool isOnlySi(false);
352  if (rhtools_.isOnlySilicon(layerId))
353  isOnlySi = true;
354 
355  for (unsigned int i = 0; i < numberOfCells; i++) {
356  bool isSi = isOnlySi || cellsOnLayer.isSi[i];
357  // initialize delta and nearest higher for i
359  float i_delta = maxDelta;
360  int i_nearestHigher = -1;
361  if (isSi) {
362  float delta = delta_c;
363  // get search box for ith hit
364  // guarantee to cover a range "outlierDeltaFactor_*delta_c"
365  auto range = outlierDeltaFactor_ * delta;
366  std::array<int, 4> search_box = lt.searchBox(
367  cellsOnLayer.x[i] - range, cellsOnLayer.x[i] + range, cellsOnLayer.y[i] - range, cellsOnLayer.y[i] + range);
368  // loop over all bins in the search box
369  for (int xBin = search_box[0]; xBin < search_box[1] + 1; ++xBin) {
370  for (int yBin = search_box[2]; yBin < search_box[3] + 1; ++yBin) {
371  // get the id of this bin
372  size_t binId = lt.getGlobalBinByBin(xBin, yBin);
373  // get the size of this bin
374  size_t binSize = lt[binId].size();
375 
376  // loop over all hits in this bin
377  for (unsigned int j = 0; j < binSize; j++) {
378  unsigned int otherId = lt[binId][j];
379  bool otherSi = isOnlySi || cellsOnLayer.isSi[otherId];
380  if (otherSi) { //silicon cells cannot talk to scintillator cells
381  float dist = distance(i, otherId, layerId, false);
382  bool foundHigher = (cellsOnLayer.rho[otherId] > cellsOnLayer.rho[i]) ||
383  (cellsOnLayer.rho[otherId] == cellsOnLayer.rho[i] &&
384  cellsOnLayer.detid[otherId] > cellsOnLayer.detid[i]);
385  // if dist == i_delta, then last comer being the nearest higher
386  if (foundHigher && dist <= i_delta) {
387  // update i_delta
388  i_delta = dist;
389  // update i_nearestHigher
390  i_nearestHigher = otherId;
391  }
392  }
393  }
394  }
395  }
396 
397  bool foundNearestHigherInSearchBox = (i_delta != maxDelta);
398  if (foundNearestHigherInSearchBox) {
399  cellsOnLayer.delta[i] = i_delta;
400  cellsOnLayer.nearestHigher[i] = i_nearestHigher;
401  } else {
402  // otherwise delta is guaranteed to be larger outlierDeltaFactor_*delta_c
403  // we can safely maximize delta to be maxDelta
404  cellsOnLayer.delta[i] = maxDelta;
405  cellsOnLayer.nearestHigher[i] = -1;
406  }
407  } else {
408  //similar to silicon
409  float delta = delta_r;
410  auto range = outlierDeltaFactor_ * delta;
411  std::array<int, 4> search_box = lt.searchBoxEtaPhi(cellsOnLayer.eta[i] - range,
412  cellsOnLayer.eta[i] + range,
413  cellsOnLayer.phi[i] - range,
414  cellsOnLayer.phi[i] + range);
415  // loop over all bins in the search box
416  for (int xBin = search_box[0]; xBin < search_box[1] + 1; ++xBin) {
417  for (int yBin = search_box[2]; yBin < search_box[3] + 1; ++yBin) {
418  // get the id of this bin
419  size_t binId = lt.getGlobalBinByBinEtaPhi(xBin, yBin);
420  // get the size of this bin
421  size_t binSize = lt[binId].size();
422 
423  // loop over all hits in this bin
424  for (unsigned int j = 0; j < binSize; j++) {
425  unsigned int otherId = lt[binId][j];
426  if (!cellsOnLayer.isSi[otherId]) { //scintillator cells cannot talk to silicon cells
427  float dist = distance(i, otherId, layerId, true);
428  bool foundHigher = (cellsOnLayer.rho[otherId] > cellsOnLayer.rho[i]) ||
429  (cellsOnLayer.rho[otherId] == cellsOnLayer.rho[i] &&
430  cellsOnLayer.detid[otherId] > cellsOnLayer.detid[i]);
431  // if dist == i_delta, then last comer being the nearest higher
432  if (foundHigher && dist <= i_delta) {
433  // update i_delta
434  i_delta = dist;
435  // update i_nearestHigher
436  i_nearestHigher = otherId;
437  }
438  }
439  }
440  }
441  }
442 
443  bool foundNearestHigherInSearchBox = (i_delta != maxDelta);
444  if (foundNearestHigherInSearchBox) {
445  cellsOnLayer.delta[i] = i_delta;
446  cellsOnLayer.nearestHigher[i] = i_nearestHigher;
447  } else {
448  // otherwise delta is guaranteed to be larger outlierDeltaFactor_*delta_c
449  // we can safely maximize delta to be maxDelta
450  cellsOnLayer.delta[i] = maxDelta;
451  cellsOnLayer.nearestHigher[i] = -1;
452  }
453  }
454  LogDebug("HGCalCLUEAlgo") << "Debugging calculateDistanceToHigher: \n"
455  << " cell: " << i << " isSilicon: " << cellsOnLayer.isSi[i]
456  << " eta: " << cellsOnLayer.eta[i] << " phi: " << cellsOnLayer.phi[i]
457  << " energy: " << cellsOnLayer.weight[i] << " density: " << cellsOnLayer.rho[i]
458  << " nearest higher: " << cellsOnLayer.nearestHigher[i]
459  << " distance: " << cellsOnLayer.delta[i] << "\n";
460  }
461 }
462 
463 template <typename T>
464 int HGCalCLUEAlgoT<T>::findAndAssignClusters(const unsigned int layerId, float delta_c, float delta_r) {
465  // this is called once per layer and endcap...
466  // so when filling the cluster temporary vector of Hexels we resize each time
467  // by the number of clusters found. This is always equal to the number of
468  // cluster centers...
469  unsigned int nClustersOnLayer = 0;
470  auto& cellsOnLayer = cells_[layerId];
471  unsigned int numberOfCells = cellsOnLayer.detid.size();
472  std::vector<int> localStack;
473  // find cluster seeds and outlier
474  for (unsigned int i = 0; i < numberOfCells; i++) {
475  float rho_c = kappa_ * cellsOnLayer.sigmaNoise[i];
476  bool isSi = rhtools_.isOnlySilicon(layerId) || cellsOnLayer.isSi[i];
477  float delta = isSi ? delta_c : delta_r;
478 
479  // initialize clusterIndex
480  cellsOnLayer.clusterIndex[i] = -1;
481  bool isSeed = (cellsOnLayer.delta[i] > delta) && (cellsOnLayer.rho[i] >= rho_c);
482  bool isOutlier = (cellsOnLayer.delta[i] > outlierDeltaFactor_ * delta) && (cellsOnLayer.rho[i] < rho_c);
483  if (isSeed) {
484  cellsOnLayer.clusterIndex[i] = nClustersOnLayer;
485  cellsOnLayer.isSeed[i] = true;
486  nClustersOnLayer++;
487  localStack.push_back(i);
488 
489  } else if (!isOutlier) {
490  cellsOnLayer.followers[cellsOnLayer.nearestHigher[i]].push_back(i);
491  }
492  }
493 
494  // need to pass clusterIndex to their followers
495  while (!localStack.empty()) {
496  int endStack = localStack.back();
497  auto& thisSeed = cellsOnLayer.followers[endStack];
498  localStack.pop_back();
499 
500  // loop over followers
501  for (int j : thisSeed) {
502  // pass id to a follower
503  cellsOnLayer.clusterIndex[j] = cellsOnLayer.clusterIndex[endStack];
504  // push this follower to localStack
505  localStack.push_back(j);
506  }
507  }
508  return nClustersOnLayer;
509 }
510 
511 template <typename T>
513  // To support the TDR geometry and also the post-TDR one (v9 onwards), we
514  // need to change the logic of the vectors containing signal to noise and
515  // thresholds. The first 3 indices will keep on addressing the different
516  // thicknesses of the Silicon detectors in CE_E , the next 3 indices will address
517  // the thicknesses of the Silicon detectors in CE_H, while the last one, number 6 (the
518  // seventh) will address the Scintillators. This change will support both
519  // geometries at the same time.
520 
521  if (initialized_)
522  return; // only need to calculate thresholds once
523 
524  initialized_ = true;
525 
526  std::vector<double> dummy;
527 
528  dummy.resize(maxNumberOfThickIndices_ + !isNose_, 0); // +1 to accomodate for the Scintillators
529  thresholds_.resize(maxlayer_, dummy);
530  v_sigmaNoise_.resize(maxlayer_, dummy);
531 
532  for (unsigned ilayer = 1; ilayer <= maxlayer_; ++ilayer) {
533  for (unsigned ithick = 0; ithick < maxNumberOfThickIndices_; ++ithick) {
534  float sigmaNoise = 0.001f * fcPerEle_ * nonAgedNoises_[ithick] * dEdXweights_[ilayer] /
535  (fcPerMip_[ithick] * thicknessCorrection_[ithick]);
536  thresholds_[ilayer - 1][ithick] = sigmaNoise * ecut_;
537  v_sigmaNoise_[ilayer - 1][ithick] = sigmaNoise;
538  LogDebug("HGCalCLUEAlgo") << "ilayer: " << ilayer << " nonAgedNoises: " << nonAgedNoises_[ithick]
539  << " fcPerEle: " << fcPerEle_ << " fcPerMip: " << fcPerMip_[ithick]
540  << " noiseMip: " << fcPerEle_ * nonAgedNoises_[ithick] / fcPerMip_[ithick]
541  << " sigmaNoise: " << sigmaNoise << "\n";
542  }
543 
544  if (!isNose_) {
545  float scintillators_sigmaNoise = 0.001f * noiseMip_ * dEdXweights_[ilayer] / sciThicknessCorrection_;
546  thresholds_[ilayer - 1][maxNumberOfThickIndices_] = ecut_ * scintillators_sigmaNoise;
547  v_sigmaNoise_[ilayer - 1][maxNumberOfThickIndices_] = scintillators_sigmaNoise;
548  LogDebug("HGCalCLUEAlgo") << "ilayer: " << ilayer << " noiseMip: " << noiseMip_
549  << " scintillators_sigmaNoise: " << scintillators_sigmaNoise << "\n";
550  }
551  }
552 }
553 
554 template <typename T>
555 void HGCalCLUEAlgoT<T>::setDensity(const unsigned int layerId) {
556  auto& cellsOnLayer = cells_[layerId];
557  unsigned int numberOfCells = cellsOnLayer.detid.size();
558  for (unsigned int i = 0; i < numberOfCells; ++i)
559  density_[cellsOnLayer.detid[i]] = cellsOnLayer.rho[i];
560 }
561 
562 template <typename T>
564  return density_;
565 }
566 
567 // explicit template instantiation
568 template class HGCalCLUEAlgoT<HGCalLayerTiles>;
569 template class HGCalCLUEAlgoT<HFNoseLayerTiles>;
FastTimerService_cff.range
range
Definition: FastTimerService_cff.py:34
mps_fire.i
i
Definition: mps_fire.py:428
PositionCalc.h
HGCScintillatorDetId::iphi
int iphi() const
get the phi index
Definition: HGCScintillatorDetId.cc:58
CaloRecHit::energy
constexpr float energy() const
Definition: CaloRecHit.h:29
photonAnalyzer_cfi.xBin
xBin
Definition: photonAnalyzer_cfi.py:81
HGCalCLUEAlgoT::findAndAssignClusters
int findAndAssignClusters(const unsigned int layerId, float delta_c, float delta_r)
Definition: HGCalCLUEAlgo.cc:464
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
etaBin
int etaBin(const l1t::HGCalMulticluster *cl)
Definition: L1EGammaEEProducer.cc:19
HGCalCLUEAlgoT::calculateLocalDensity
void calculateLocalDensity(const TILE &lt, const unsigned int layerId, float delta_c, float delta_r)
Definition: HGCalCLUEAlgo.cc:247
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
HGCScintillatorDetId::ieta
int ieta() const
Definition: HGCScintillatorDetId.h:57
reco::CaloID::DET_HGCAL_ENDCAP
Definition: CaloID.h:30
BeamMonitor_cff.phiBin
phiBin
Definition: BeamMonitor_cff.py:77
hitfit::delta_r
double delta_r(const Fourvec &a, const Fourvec &b)
Find the distance between two four-vectors in the two-dimensional space .
Definition: fourvec.cc:238
DetId::det
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
reco::BasicCluster
CaloCluster BasicCluster
Definition: BasicClusterFwd.h:13
HGCalCLUEAlgoT
Definition: HGCalCLUEAlgo.h:30
HGCalCLUEAlgoT::getEventSetupPerAlgorithm
void getEventSetupPerAlgorithm(const edm::EventSetup &es) override
Definition: HGCalCLUEAlgo.cc:19
CaloID.h
edm::SortedCollection
Definition: SortedCollection.h:49
photonAnalyzer_cfi.yBin
yBin
Definition: photonAnalyzer_cfi.py:85
findQualityFiles.v
v
Definition: findQualityFiles.py:179
CaloRecHit::detid
constexpr const DetId & detid() const
Definition: CaloRecHit.h:33
allConversions_cfi.maxDelta
maxDelta
Definition: allConversions_cfi.py:44
python.cmstools.all
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
GetRecoTauVFromDQM_MC_cff.cl
cl
Definition: GetRecoTauVFromDQM_MC_cff.py:38
DetId
Definition: DetId.h:17
DetId::HGCalHSi
Definition: DetId.h:33
HGCalCLUEAlgoT::makeClusters
void makeClusters() override
Definition: HGCalCLUEAlgo.cc:99
HGCalCLUEAlgoT::getClusters
std::vector< reco::BasicCluster > getClusters(bool) override
Definition: HGCalCLUEAlgo.cc:131
HGCalCLUEAlgoT::populate
void populate(const HGCRecHitCollection &hits) override
Definition: HGCalCLUEAlgo.cc:27
hgcal_clustering::Density
std::map< DetId, float > Density
Definition: HGCalClusteringAlgoBase.h:44
PVValHelper::eta
Definition: PVValidationHelpers.h:70
HGCalCLUEAlgoT::prepareDataStructures
void prepareDataStructures(const unsigned int layerId)
Definition: HGCalCLUEAlgo.cc:79
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
hgcal_clustering
Definition: HGCalClusteringAlgoBase.h:17
Point3DBase< float, GlobalTag >
HGCalCLUEAlgoT::computeThreshold
void computeThreshold()
Definition: HGCalCLUEAlgo.cc:512
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
HGCRecHit
Definition: HGCRecHit.h:14
HGCalCLUEAlgoT::getDensity
Density getDensity() override
Definition: HGCalCLUEAlgo.cc:563
CaloSubdetectorGeometry.h
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
dumpMFGeometry_cfg.delta
delta
Definition: dumpMFGeometry_cfg.py:25
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
IdealGeometryRecord.h
HGCalCLUEAlgoT::setDensity
void setDensity(const unsigned int layerId)
Definition: HGCalCLUEAlgo.cc:555
edm::EventSetup
Definition: EventSetup.h:58
HcalSubdetector.h
cmsLHEtoEOSManager.l
l
Definition: cmsLHEtoEOSManager.py:204
CaloCellGeometry.h
HGCalCLUEAlgoT::calculatePosition
math::XYZPoint calculatePosition(const std::vector< int > &v, const unsigned int layerId) const
Definition: HGCalCLUEAlgo.cc:187
AlignmentPI::calculatePosition
std::pair< double, double > calculatePosition(TVirtualPad *myPad, int boundary)
Definition: AlignmentPayloadInspectorHelper.h:770
HGCalCLUEAlgo.h
HGCScintillatorDetId
Definition: HGCScintillatorDetId.h:23
T
long double T
Definition: Basic3DVectorLD.h:48
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
unpackBuffers-CaloStage1.offsets
offsets
Definition: unpackBuffers-CaloStage1.py:127
L1TowerCalibrationProducer_cfi.iEta
iEta
Definition: L1TowerCalibrationProducer_cfi.py:60
HGCHEF
Definition: ForwardSubdetector.h:9
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
dummy
Definition: DummySelector.h:38
HGCalCLUEAlgoT::calculateDistanceToHigher
void calculateDistanceToHigher(const TILE &lt, const unsigned int layerId, float delta_c, float delta_r)
Definition: HGCalCLUEAlgo.cc:345
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
hltrates_dqm_sourceclient-live_cfg.offset
offset
Definition: hltrates_dqm_sourceclient-live_cfg.py:82
HLT_FULL_cff.distance
distance
Definition: HLT_FULL_cff.py:7746