CMS 3D CMS Logo

HGCalCLUEAlgo.cc
Go to the documentation of this file.
2 
3 // Geometry
8 
10 //
12 #include "oneapi/tbb/task_arena.h"
13 #include "oneapi/tbb.h"
14 #include <limits>
16 
17 using namespace hgcal_clustering;
18 
19 template <typename T, typename STRATEGY>
21  cells_.clear();
22  numberOfClustersPerLayer_.clear();
23  cells_.resize(2 * (maxlayer_ + 1));
24  numberOfClustersPerLayer_.resize(2 * (maxlayer_ + 1), 0);
25 }
26 
27 template <typename T, typename STRATEGY>
29  // loop over all hits and create the Hexel structure, skip energies below ecut
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  if constexpr (std::is_same_v<STRATEGY, HGCalScintillatorStrategy>) {
67  cells_[layer].dim1.emplace_back(position.eta());
68  cells_[layer].dim2.emplace_back(position.phi());
69  } // else, isSilicon == true and eta phi values will not be used
70  else {
71  cells_[layer].dim1.emplace_back(position.x());
72  cells_[layer].dim2.emplace_back(position.y());
73  }
74  cells_[layer].weight.emplace_back(hgrh.energy());
75  cells_[layer].sigmaNoise.emplace_back(sigmaNoise);
76  }
77 }
78 
79 template <typename T, typename STRATEGY>
81  auto cellsSize = cells_[l].detid.size();
82  cells_[l].rho.resize(cellsSize, 0.f);
83  cells_[l].delta.resize(cellsSize, 9999999);
84  cells_[l].nearestHigher.resize(cellsSize, -1);
85  cells_[l].clusterIndex.resize(cellsSize, -1);
86  cells_[l].followers.resize(cellsSize);
87  cells_[l].isSeed.resize(cellsSize, false);
88 }
89 
90 // Create a vector of Hexels associated to one cluster from a collection of
91 // HGCalRecHits - this can be used directly to make the final cluster list -
92 // this method can be invoked multiple times for the same event with different
93 // input (reset should be called between events)
94 template <typename T, typename STRATEGY>
96  // assign all hits in each layer to a cluster core
97  tbb::this_task_arena::isolate([&] {
98  tbb::parallel_for(size_t(0), size_t(2 * maxlayer_ + 2), [&](size_t i) {
99  prepareDataStructures(i);
100  T lt;
101  lt.clear();
102  lt.fill(cells_[i].dim1, cells_[i].dim2);
103 
104  float delta;
105  if constexpr (std::is_same_v<STRATEGY, HGCalSiliconStrategy>) {
106  // maximum search distance (critical distance) for local density calculation
107  float delta_c;
108  if (i % maxlayer_ < lastLayerEE_)
109  delta_c = vecDeltas_[0];
110  else if (i % maxlayer_ < (firstLayerBH_ - 1))
111  delta_c = vecDeltas_[1];
112  else
113  delta_c = vecDeltas_[2];
114  delta = delta_c;
115  } else {
116  float delta_r = vecDeltas_[3];
117  delta = delta_r;
118  }
119  LogDebug("HGCalCLUEAlgo") << "maxlayer: " << maxlayer_ << " lastLayerEE: " << lastLayerEE_
120  << " firstLayerBH: " << firstLayerBH_ << "\n";
121 
122  calculateLocalDensity(lt, i, delta);
123  calculateDistanceToHigher(lt, i, delta);
124  numberOfClustersPerLayer_[i] = findAndAssignClusters(i, delta);
125  });
126  });
127 }
128 
129 template <typename T, typename STRATEGY>
130 std::vector<reco::BasicCluster> HGCalCLUEAlgoT<T, STRATEGY>::getClusters(bool) {
131  std::vector<int> offsets(numberOfClustersPerLayer_.size(), 0);
132 
133  int maxClustersOnLayer = numberOfClustersPerLayer_[0];
134 
135  for (unsigned layerId = 1; layerId < offsets.size(); ++layerId) {
136  offsets[layerId] = offsets[layerId - 1] + numberOfClustersPerLayer_[layerId - 1];
137 
138  maxClustersOnLayer = std::max(maxClustersOnLayer, numberOfClustersPerLayer_[layerId]);
139  }
140 
141  auto totalNumberOfClusters = offsets.back() + numberOfClustersPerLayer_.back();
142  clusters_v_.resize(totalNumberOfClusters);
143  std::vector<std::vector<int>> cellsIdInCluster;
144  cellsIdInCluster.reserve(maxClustersOnLayer);
145 
146  for (unsigned int layerId = 0; layerId < 2 * maxlayer_ + 2; ++layerId) {
147  cellsIdInCluster.resize(numberOfClustersPerLayer_[layerId]);
148  auto& cellsOnLayer = cells_[layerId];
149  unsigned int numberOfCells = cellsOnLayer.detid.size();
150  auto firstClusterIdx = offsets[layerId];
151 
152  for (unsigned int i = 0; i < numberOfCells; ++i) {
153  auto clusterIndex = cellsOnLayer.clusterIndex[i];
154  if (clusterIndex != -1)
155  cellsIdInCluster[clusterIndex].push_back(i);
156  }
157 
158  std::vector<std::pair<DetId, float>> thisCluster;
159 
160  for (auto& cl : cellsIdInCluster) {
162  float energy = 0.f;
163  int seedDetId = -1;
164 
165  for (auto cellIdx : cl) {
166  energy += cellsOnLayer.weight[cellIdx];
167  thisCluster.emplace_back(cellsOnLayer.detid[cellIdx], 1.f);
168  if (cellsOnLayer.isSeed[cellIdx]) {
169  seedDetId = cellsOnLayer.detid[cellIdx];
170  }
171  }
172  auto globalClusterIndex = cellsOnLayer.clusterIndex[cl[0]] + firstClusterIdx;
173 
174  clusters_v_[globalClusterIndex] =
176  clusters_v_[globalClusterIndex].setSeed(seedDetId);
177  thisCluster.clear();
178  }
179 
180  cellsIdInCluster.clear();
181  }
182  return clusters_v_;
183 }
184 template <typename T, typename STRATEGY>
186  const unsigned int layerId,
187  float delta,
189  auto& cellsOnLayer = cells_[layerId];
190  unsigned int numberOfCells = cellsOnLayer.detid.size();
191  for (unsigned int i = 0; i < numberOfCells; i++) {
192  std::array<int, 4> search_box = lt.searchBox(cellsOnLayer.dim1[i] - delta,
193  cellsOnLayer.dim1[i] + delta,
194  cellsOnLayer.dim2[i] - delta,
195  cellsOnLayer.dim2[i] + delta);
196 
197  for (int xBin = search_box[0]; xBin < search_box[1] + 1; ++xBin) {
198  for (int yBin = search_box[2]; yBin < search_box[3] + 1; ++yBin) {
199  int binId = lt.getGlobalBinByBin(xBin, yBin);
200  size_t binSize = lt[binId].size();
201 
202  for (unsigned int j = 0; j < binSize; j++) {
203  unsigned int otherId = lt[binId][j];
204  if (distance(lt, i, otherId, layerId) < delta) {
205  cellsOnLayer.rho[i] += (i == otherId ? 1.f : 0.5f) * cellsOnLayer.weight[otherId];
206  }
207  }
208  }
209  }
210  LogDebug("HGCalCLUEAlgo") << "Debugging calculateLocalDensity: \n"
211  << " cell: " << i << " eta: " << cellsOnLayer.dim1[i] << " phi: " << cellsOnLayer.dim2[i]
212  << " energy: " << cellsOnLayer.weight[i] << " density: " << cellsOnLayer.rho[i] << "\n";
213  }
214 }
215 template <typename T, typename STRATEGY>
217  const unsigned int layerId,
218  float delta,
220  auto& cellsOnLayer = cells_[layerId];
221  unsigned int numberOfCells = cellsOnLayer.detid.size();
222  for (unsigned int i = 0; i < numberOfCells; i++) {
223  std::array<int, 4> search_box = lt.searchBox(cellsOnLayer.dim1[i] - delta,
224  cellsOnLayer.dim1[i] + delta,
225  cellsOnLayer.dim2[i] - delta,
226  cellsOnLayer.dim2[i] + delta);
227  cellsOnLayer.rho[i] += cellsOnLayer.weight[i];
228  float northeast(0), northwest(0), southeast(0), southwest(0), all(0);
229  for (int etaBin = search_box[0]; etaBin < search_box[1] + 1; ++etaBin) {
230  for (int phiBin = search_box[2]; phiBin < search_box[3] + 1; ++phiBin) {
231  int phi = (phiBin % T::type::nRows);
232  int binId = lt.getGlobalBinByBin(etaBin, phi);
233  size_t binSize = lt[binId].size();
234  for (unsigned int j = 0; j < binSize; j++) {
235  unsigned int otherId = lt[binId][j];
236  if (distance(lt, i, otherId, layerId) < delta) {
237  int iPhi = HGCScintillatorDetId(cellsOnLayer.detid[i]).iphi();
238  int otherIPhi = HGCScintillatorDetId(cellsOnLayer.detid[otherId]).iphi();
239  int iEta = HGCScintillatorDetId(cellsOnLayer.detid[i]).ieta();
240  int otherIEta = HGCScintillatorDetId(cellsOnLayer.detid[otherId]).ieta();
241  int dIPhi = otherIPhi - iPhi;
242  dIPhi += abs(dIPhi) < 2 ? 0
243  : dIPhi < 0 ? scintMaxIphi_
244  : -scintMaxIphi_; // cells with iPhi=288 and iPhi=1 should be neiboring cells
245  int dIEta = otherIEta - iEta;
246  LogDebug("HGCalCLUEAlgo") << " Debugging calculateLocalDensity for Scintillator: \n"
247  << " cell: " << otherId << " energy: " << cellsOnLayer.weight[otherId]
248  << " otherIPhi: " << otherIPhi << " iPhi: " << iPhi << " otherIEta: " << otherIEta
249  << " iEta: " << iEta << "\n";
250 
251  if (otherId != i) {
252  auto neighborCellContribution = 0.5f * cellsOnLayer.weight[otherId];
253  all += neighborCellContribution;
254  if (dIPhi >= 0 && dIEta >= 0)
255  northeast += neighborCellContribution;
256  if (dIPhi <= 0 && dIEta >= 0)
257  southeast += neighborCellContribution;
258  if (dIPhi >= 0 && dIEta <= 0)
259  northwest += neighborCellContribution;
260  if (dIPhi <= 0 && dIEta <= 0)
261  southwest += neighborCellContribution;
262  }
263  LogDebug("HGCalCLUEAlgo") << " Debugging calculateLocalDensity for Scintillator: \n"
264  << " northeast: " << northeast << " southeast: " << southeast
265  << " northwest: " << northwest << " southwest: " << southwest << "\n";
266  }
267  }
268  }
269  }
270  float neighborsval = (std::max(northeast, northwest) > std::max(southeast, southwest))
271  ? std::max(northeast, northwest)
272  : std::max(southeast, southwest);
273  if (use2x2_)
274  cellsOnLayer.rho[i] += neighborsval;
275  else
276  cellsOnLayer.rho[i] += all;
277  LogDebug("HGCalCLUEAlgo") << "Debugging calculateLocalDensity: \n"
278  << " cell: " << i << " eta: " << cellsOnLayer.dim1[i] << " phi: " << cellsOnLayer.dim2[i]
279  << " energy: " << cellsOnLayer.weight[i] << " density: " << cellsOnLayer.rho[i] << "\n";
280  }
281 }
282 template <typename T, typename STRATEGY>
283 void HGCalCLUEAlgoT<T, STRATEGY>::calculateLocalDensity(const T& lt, const unsigned int layerId, float delta) {
284  if constexpr (std::is_same_v<STRATEGY, HGCalSiliconStrategy>) {
285  calculateLocalDensity(lt, layerId, delta, HGCalSiliconStrategy());
286  } else {
287  calculateLocalDensity(lt, layerId, delta, HGCalScintillatorStrategy());
288  }
289 }
290 
291 template <typename T, typename STRATEGY>
292 void HGCalCLUEAlgoT<T, STRATEGY>::calculateDistanceToHigher(const T& lt, const unsigned int layerId, float delta) {
293  auto& cellsOnLayer = cells_[layerId];
294  unsigned int numberOfCells = cellsOnLayer.detid.size();
295 
296  for (unsigned int i = 0; i < numberOfCells; i++) {
297  // initialize delta and nearest higher for i
299  float i_delta = maxDelta;
300  int i_nearestHigher = -1;
301  auto range = outlierDeltaFactor_ * delta;
302  std::array<int, 4> search_box = lt.searchBox(cellsOnLayer.dim1[i] - range,
303  cellsOnLayer.dim1[i] + range,
304  cellsOnLayer.dim2[i] - range,
305  cellsOnLayer.dim2[i] + range);
306  // loop over all bins in the search box
307  for (int dim1Bin = search_box[0]; dim1Bin < search_box[1] + 1; ++dim1Bin) {
308  for (int dim2Bin = search_box[2]; dim2Bin < search_box[3] + 1; ++dim2Bin) {
309  // get the id of this bin
310  size_t binId = lt.getGlobalBinByBin(dim1Bin, dim2Bin);
311  if constexpr (std::is_same_v<STRATEGY, HGCalScintillatorStrategy>)
312  binId = lt.getGlobalBinByBin(dim1Bin, (dim2Bin % T::type::nRows));
313  // get the size of this bin
314  size_t binSize = lt[binId].size();
315 
316  // loop over all hits in this bin
317  for (unsigned int j = 0; j < binSize; j++) {
318  unsigned int otherId = lt[binId][j];
319  float dist = distance(lt, i, otherId, layerId);
320  bool foundHigher =
321  (cellsOnLayer.rho[otherId] > cellsOnLayer.rho[i]) ||
322  (cellsOnLayer.rho[otherId] == cellsOnLayer.rho[i] && cellsOnLayer.detid[otherId] > cellsOnLayer.detid[i]);
323  if (foundHigher && dist <= i_delta) {
324  // update i_delta
325  i_delta = dist;
326  // update i_nearestHigher
327  i_nearestHigher = otherId;
328  }
329  }
330  }
331  }
332  bool foundNearestHigherInSearchBox = (i_delta != maxDelta);
333  if (foundNearestHigherInSearchBox) {
334  cellsOnLayer.delta[i] = i_delta;
335  cellsOnLayer.nearestHigher[i] = i_nearestHigher;
336  } else {
337  // otherwise delta is guaranteed to be larger outlierDeltaFactor_*delta_c
338  // we can safely maximize delta to be maxDelta
339  cellsOnLayer.delta[i] = maxDelta;
340  cellsOnLayer.nearestHigher[i] = -1;
341  }
342 
343  LogDebug("HGCalCLUEAlgo") << "Debugging calculateDistanceToHigher: \n"
344  << " cell: " << i << " eta: " << cellsOnLayer.dim1[i] << " phi: " << cellsOnLayer.dim2[i]
345  << " energy: " << cellsOnLayer.weight[i] << " density: " << cellsOnLayer.rho[i]
346  << " nearest higher: " << cellsOnLayer.nearestHigher[i]
347  << " distance: " << cellsOnLayer.delta[i] << "\n";
348  }
349 }
350 
351 template <typename T, typename STRATEGY>
352 int HGCalCLUEAlgoT<T, STRATEGY>::findAndAssignClusters(const unsigned int layerId, float delta) {
353  // this is called once per layer and endcap...
354  // so when filling the cluster temporary vector of Hexels we resize each time
355  // by the number of clusters found. This is always equal to the number of
356  // cluster centers...
357  unsigned int nClustersOnLayer = 0;
358  auto& cellsOnLayer = cells_[layerId];
359  unsigned int numberOfCells = cellsOnLayer.detid.size();
360  std::vector<int> localStack;
361  // find cluster seeds and outlier
362  for (unsigned int i = 0; i < numberOfCells; i++) {
363  float rho_c = kappa_ * cellsOnLayer.sigmaNoise[i];
364  // initialize clusterIndex
365  cellsOnLayer.clusterIndex[i] = -1;
366  bool isSeed = (cellsOnLayer.delta[i] > delta) && (cellsOnLayer.rho[i] >= rho_c);
367  bool isOutlier = (cellsOnLayer.delta[i] > outlierDeltaFactor_ * delta) && (cellsOnLayer.rho[i] < rho_c);
368  if (isSeed) {
369  cellsOnLayer.clusterIndex[i] = nClustersOnLayer;
370  cellsOnLayer.isSeed[i] = true;
371  nClustersOnLayer++;
372  localStack.push_back(i);
373 
374  } else if (!isOutlier) {
375  cellsOnLayer.followers[cellsOnLayer.nearestHigher[i]].push_back(i);
376  }
377  }
378 
379  // need to pass clusterIndex to their followers
380  while (!localStack.empty()) {
381  int endStack = localStack.back();
382  auto& thisSeed = cellsOnLayer.followers[endStack];
383  localStack.pop_back();
384 
385  // loop over followers
386  for (int j : thisSeed) {
387  // pass id to a follower
388  cellsOnLayer.clusterIndex[j] = cellsOnLayer.clusterIndex[endStack];
389  // push this follower to localStack
390  localStack.push_back(j);
391  }
392  }
393  return nClustersOnLayer;
394 }
395 
396 template <typename T, typename STRATEGY>
398  // To support the TDR geometry and also the post-TDR one (v9 onwards), we
399  // need to change the logic of the vectors containing signal to noise and
400  // thresholds. The first 3 indices will keep on addressing the different
401  // thicknesses of the Silicon detectors in CE_E , the next 3 indices will address
402  // the thicknesses of the Silicon detectors in CE_H, while the last one, number 6 (the
403  // seventh) will address the Scintillators. This change will support both
404  // geometries at the same time.
405 
406  if (initialized_)
407  return; // only need to calculate thresholds once
408 
409  initialized_ = true;
410 
411  std::vector<double> dummy;
412 
413  dummy.resize(maxNumberOfThickIndices_ + !isNose_, 0); // +1 to accomodate for the Scintillators
414  thresholds_.resize(maxlayer_, dummy);
415  v_sigmaNoise_.resize(maxlayer_, dummy);
416 
417  for (unsigned ilayer = 1; ilayer <= maxlayer_; ++ilayer) {
418  for (unsigned ithick = 0; ithick < maxNumberOfThickIndices_; ++ithick) {
419  float sigmaNoise = 0.001f * fcPerEle_ * nonAgedNoises_[ithick] * dEdXweights_[ilayer] /
420  (fcPerMip_[ithick] * thicknessCorrection_[ithick]);
421  thresholds_[ilayer - 1][ithick] = sigmaNoise * ecut_;
422  v_sigmaNoise_[ilayer - 1][ithick] = sigmaNoise;
423  LogDebug("HGCalCLUEAlgo") << "ilayer: " << ilayer << " nonAgedNoises: " << nonAgedNoises_[ithick]
424  << " fcPerEle: " << fcPerEle_ << " fcPerMip: " << fcPerMip_[ithick]
425  << " noiseMip: " << fcPerEle_ * nonAgedNoises_[ithick] / fcPerMip_[ithick]
426  << " sigmaNoise: " << sigmaNoise << "\n";
427  }
428 
429  if (!isNose_) {
430  float scintillators_sigmaNoise = 0.001f * noiseMip_ * dEdXweights_[ilayer] / sciThicknessCorrection_;
431  thresholds_[ilayer - 1][maxNumberOfThickIndices_] = ecut_ * scintillators_sigmaNoise;
432  v_sigmaNoise_[ilayer - 1][maxNumberOfThickIndices_] = scintillators_sigmaNoise;
433  LogDebug("HGCalCLUEAlgo") << "ilayer: " << ilayer << " noiseMip: " << noiseMip_
434  << " scintillators_sigmaNoise: " << scintillators_sigmaNoise << "\n";
435  }
436  }
437 }
438 
439 // explicit template instantiation
constexpr int ieta() const
void computeThreshold()
constexpr int iphi() const
get the phi index
constexpr const DetId & detid() const
Definition: CaloRecHit.h:33
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
CaloCluster BasicCluster
void calculateDistanceToHigher(const TILE &lt, const unsigned int layerId, float delta)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
void populate(const HGCRecHitCollection &hits) override
constexpr float energy() const
Definition: CaloRecHit.h:29
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void makeClusters() override
double f[11][100]
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Definition: DetId.h:17
void getEventSetupPerAlgorithm(const edm::EventSetup &es) override
void calculateLocalDensity(const TILE &lt, const unsigned int layerId, float delta)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void prepareDataStructures(const unsigned int layerId)
static int position[264][3]
Definition: ReadPGInfo.cc:289
long double T
strategy
Definition: nnet_common.h:18
int findAndAssignClusters(const unsigned int layerId, float delta)
std::vector< reco::BasicCluster > getClusters(bool) override
#define LogDebug(id)