CMS 3D CMS Logo

HGCalImagingAlgo.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 
16  // loop over all hits and create the Hexel structure, skip energies below ecut
17 
18  if (dependSensor_) {
19  // for each layer and wafer calculate the thresholds (sigmaNoise and energy)
20  // once
22  }
23 
24  std::vector<bool> firstHit(2 * (maxlayer + 1), true);
25  for (unsigned int i = 0; i < hits.size(); ++i) {
26 
27  const HGCRecHit &hgrh = hits[i];
28  DetId detid = hgrh.detid();
29  unsigned int layer = rhtools_.getLayerWithOffset(detid);
30  float thickness = 0.f;
31  // set sigmaNoise default value 1 to use kappa value directly in case of
32  // sensor-independent thresholds
33  float sigmaNoise = 1.f;
34  if (dependSensor_) {
35  thickness = rhtools_.getSiThickness(detid);
36  int thickness_index = rhtools_.getSiThickIndex(detid);
37  if (thickness_index == -1)
38  thickness_index = 3;
39  double storedThreshold = thresholds_[layer - 1][thickness_index];
40  sigmaNoise = v_sigmaNoise_[layer - 1][thickness_index];
41 
42  if (hgrh.energy() < storedThreshold)
43  continue; // this sets the ZS threshold at ecut times the sigma noise
44  // for the sensor
45  }
46  if (!dependSensor_ && hgrh.energy() < ecut_)
47  continue;
48 
49  // map layers from positive endcap (z) to layer + maxlayer+1 to prevent
50  // mixing up hits from different sides
51  layer += int(rhtools_.zside(detid) > 0) * (maxlayer + 1);
52 
53  // determine whether this is a half-hexagon
54  bool isHalf = rhtools_.isHalfCell(detid);
56 
57  // here's were the KDNode is passed its dims arguments - note that these are
58  // *copied* from the Hexel
59  points_[layer].emplace_back(
60  Hexel(hgrh, detid, isHalf, sigmaNoise, thickness, &rhtools_),
61  position.x(), position.y());
62 
63  // for each layer, store the minimum and maximum x and y coordinates for the
64  // KDTreeBox boundaries
65  if (firstHit[layer]) {
66  minpos_[layer][0] = position.x();
67  minpos_[layer][1] = position.y();
68  maxpos_[layer][0] = position.x();
69  maxpos_[layer][1] = position.y();
70  firstHit[layer] = false;
71  } else {
72  minpos_[layer][0] = std::min((float)position.x(), minpos_[layer][0]);
73  minpos_[layer][1] = std::min((float)position.y(), minpos_[layer][1]);
74  maxpos_[layer][0] = std::max((float)position.x(), maxpos_[layer][0]);
75  maxpos_[layer][1] = std::max((float)position.y(), maxpos_[layer][1]);
76  }
77 
78  } // end loop hits
79 }
80 // Create a vector of Hexels associated to one cluster from a collection of
81 // HGCalRecHits - this can be used directly to make the final cluster list -
82 // this method can be invoked multiple times for the same event with different
83 // input (reset should be called between events)
85  layerClustersPerLayer_.resize(2 * maxlayer + 2);
86  // assign all hits in each layer to a cluster core or halo
87  tbb::this_task_arena::isolate([&] {
88  tbb::parallel_for(size_t(0), size_t(2 * maxlayer + 2), [&](size_t i) {
89  KDTreeBox bounds(minpos_[i][0], maxpos_[i][0], minpos_[i][1], maxpos_[i][1]);
90  KDTree hit_kdtree;
91  hit_kdtree.build(points_[i], bounds);
92 
93  unsigned int actualLayer =
94  i > maxlayer
95  ? (i - (maxlayer + 1))
96  : i; // maps back from index used for KD trees to actual layer
97 
98  double maxdensity = calculateLocalDensity(
99  points_[i], hit_kdtree, actualLayer); // also stores rho (energy
100  // density) for each point (node)
101  // calculate distance to nearest point with higher density storing
102  // distance (delta) and point's index
104  findAndAssignClusters(points_[i], hit_kdtree, maxdensity, bounds,
105  actualLayer, layerClustersPerLayer_[i]);
106  });
107  });
108 }
109 
110 std::vector<reco::BasicCluster> HGCalImagingAlgo::getClusters(bool doSharing) {
111 
113  std::vector<std::pair<DetId, float>> thisCluster;
114  for (auto &clsOnLayer : layerClustersPerLayer_) {
115  for (unsigned int i = 0; i < clsOnLayer.size(); ++i) {
116  double energy = 0;
117  Point position;
118 
119  if (doSharing) {
120 
121  std::vector<unsigned> seeds = findLocalMaximaInCluster(clsOnLayer[i]);
122  // sharing found seeds.size() sub-cluster seeds in cluster i
123 
124  std::vector<std::vector<double>> fractions;
125  // first pass can have noise it in
126  shareEnergy(clsOnLayer[i], seeds, fractions);
127 
128  // reset and run second pass after vetoing seeds
129  // that result in trivial clusters (less than 2 effective cells)
130 
131  for (unsigned isub = 0; isub < fractions.size(); ++isub) {
132  double effective_hits = 0.0;
133  double energy =
134  calculateEnergyWithFraction(clsOnLayer[i], fractions[isub]);
135  Point position =
136  calculatePositionWithFraction(clsOnLayer[i], fractions[isub]);
137 
138  for (unsigned ihit = 0; ihit < fractions[isub].size(); ++ihit) {
139  const double fraction = fractions[isub][ihit];
140  if (fraction > 1e-7) {
141  effective_hits += fraction;
142  thisCluster.emplace_back(clsOnLayer[i][ihit].data.detid,
143  fraction);
144  }
145  }
146 
147  if (verbosity_ < pINFO) {
148  std::cout << "\t******** NEW CLUSTER (SHARING) ********"
149  << std::endl;
150  std::cout << "\tEff. No. of cells = " << effective_hits
151  << std::endl;
152  std::cout << "\t Energy = " << energy << std::endl;
153  std::cout << "\t Phi = " << position.phi()
154  << std::endl;
155  std::cout << "\t Eta = " << position.eta()
156  << std::endl;
157  std::cout << "\t*****************************" << std::endl;
158  }
159  clusters_v_.emplace_back(energy, position, caloID, thisCluster,
160  algoId_);
161  thisCluster.clear();
162  }
163  } else {
164  position = calculatePosition(clsOnLayer[i]); // energy-weighted position
165  // std::vector< KDNode >::iterator it;
166  for (auto &it : clsOnLayer[i]) {
167  energy += it.data.isHalo ? 0. : it.data.weight;
168  // use fraction to store whether this is a Halo hit or not
169  thisCluster.emplace_back(it.data.detid, (it.data.isHalo ? 0.f : 1.f));
170  }
171  if (verbosity_ < pINFO) {
172  std::cout << "******** NEW CLUSTER (HGCIA) ********" << std::endl;
173  std::cout << "Index " << i << std::endl;
174  std::cout << "No. of cells = " << clsOnLayer[i].size() << std::endl;
175  std::cout << " Energy = " << energy << std::endl;
176  std::cout << " Phi = " << position.phi() << std::endl;
177  std::cout << " Eta = " << position.eta() << std::endl;
178  std::cout << "*****************************" << std::endl;
179  }
180  clusters_v_.emplace_back(energy, position, caloID, thisCluster, algoId_);
181  thisCluster.clear();
182  }
183  }
184  }
185  return clusters_v_;
186 }
187 
189 HGCalImagingAlgo::calculatePosition(std::vector<KDNode> &v) const {
190  float total_weight = 0.f;
191  float x = 0.f;
192  float y = 0.f;
193  float z = 0.f;
194  unsigned int v_size = v.size();
195  unsigned int maxEnergyIndex = 0;
196  float maxEnergyValue = 0;
197  bool haloOnlyCluster = true;
198 
199  // loop over hits in cluster candidate building up weight for
200  // energy-weighted position calculation and determining the maximum
201  // energy hit in case this is a halo-only cluster
202  for (unsigned int i = 0; i < v_size; i++) {
203  if (!v[i].data.isHalo) {
204  haloOnlyCluster = false;
205  total_weight += v[i].data.weight;
206  x += v[i].data.x * v[i].data.weight;
207  y += v[i].data.y * v[i].data.weight;
208  z += v[i].data.z * v[i].data.weight;
209  } else {
210  if (v[i].data.weight > maxEnergyValue) {
211  maxEnergyValue = v[i].data.weight;
212  maxEnergyIndex = i;
213  }
214  }
215  }
216 
217  if (!haloOnlyCluster) {
218  if (total_weight != 0) {
219  auto inv_tot_weight = 1. / total_weight;
220  return math::XYZPoint(x * inv_tot_weight, y * inv_tot_weight,
221  z * inv_tot_weight);
222  }
223  } else if (v_size > 0) {
224  // return position of hit with maximum energy
225  return math::XYZPoint(v[maxEnergyIndex].data.x, v[maxEnergyIndex].data.y,
226  v[maxEnergyIndex].data.z);
227  }
228  return math::XYZPoint(0, 0, 0);
229 }
230 
231 double HGCalImagingAlgo::calculateLocalDensity(std::vector<KDNode> &nd,
232  KDTree &lp,
233  const unsigned int layer) const {
234 
235  double maxdensity = 0.;
236  float delta_c; // maximum search distance (critical distance) for local
237  // density calculation
238  if (layer <= lastLayerEE)
239  delta_c = vecDeltas_[0];
240  else if (layer <= lastLayerFH)
241  delta_c = vecDeltas_[1];
242  else
243  delta_c = vecDeltas_[2];
244 
245  // for each node calculate local density rho and store it
246  for (unsigned int i = 0; i < nd.size(); ++i) {
247  // speec up search by looking within +/- delta_c window only
248  KDTreeBox search_box(nd[i].dims[0] - delta_c, nd[i].dims[0] + delta_c,
249  nd[i].dims[1] - delta_c, nd[i].dims[1] + delta_c);
250  std::vector<KDNode> found;
251  lp.search(search_box, found);
252  const unsigned int found_size = found.size();
253  for (unsigned int j = 0; j < found_size; j++) {
254  if (distance(nd[i].data, found[j].data) < delta_c) {
255  nd[i].data.rho += found[j].data.weight;
256  maxdensity = std::max(maxdensity, nd[i].data.rho);
257  }
258  } // end loop found
259  } // end loop nodes
260  return maxdensity;
261 }
262 
263 double
264 HGCalImagingAlgo::calculateDistanceToHigher(std::vector<KDNode> &nd) const {
265 
266  // sort vector of Hexels by decreasing local density
267  std::vector<size_t> rs = sorted_indices(nd);
268 
269  double maxdensity = 0.0;
270  int nearestHigher = -1;
271 
272  if (!rs.empty())
273  maxdensity = nd[rs[0]].data.rho;
274  else
275  return maxdensity; // there are no hits
276  double dist2 = 0.;
277  // start by setting delta for the highest density hit to
278  // the most distant hit - this is a convention
279 
280  for (auto &j : nd) {
281  double tmp = distance2(nd[rs[0]].data, j.data);
282  if (tmp > dist2)
283  dist2 = tmp;
284  }
285  nd[rs[0]].data.delta = std::sqrt(dist2);
286  nd[rs[0]].data.nearestHigher = nearestHigher;
287 
288  // now we save the largest distance as a starting point
289  const double max_dist2 = dist2;
290  const unsigned int nd_size = nd.size();
291 
292  for (unsigned int oi = 1; oi < nd_size;
293  ++oi) { // start from second-highest density
294  dist2 = max_dist2;
295  unsigned int i = rs[oi];
296  // we only need to check up to oi since hits
297  // are ordered by decreasing density
298  // and all points coming BEFORE oi are guaranteed to have higher rho
299  // and the ones AFTER to have lower rho
300  for (unsigned int oj = 0; oj < oi; ++oj) {
301  unsigned int j = rs[oj];
302  double tmp = distance2(nd[i].data, nd[j].data);
303  if (tmp <= dist2) { // this "<=" instead of "<" addresses the (rare) case
304  // when there are only two hits
305  dist2 = tmp;
306  nearestHigher = j;
307  }
308  }
309  nd[i].data.delta = std::sqrt(dist2);
310  nd[i].data.nearestHigher =
311  nearestHigher; // this uses the original unsorted hitlist
312  }
313  return maxdensity;
314 }
316  std::vector<KDNode> &nd, KDTree &lp, double maxdensity, KDTreeBox &bounds,
317  const unsigned int layer,
318  std::vector<std::vector<KDNode>> &clustersOnLayer) const {
319 
320  // this is called once per layer and endcap...
321  // so when filling the cluster temporary vector of Hexels we resize each time
322  // by the number of clusters found. This is always equal to the number of
323  // cluster centers...
324 
325  unsigned int nClustersOnLayer = 0;
326  float delta_c; // critical distance
327  if (layer <= lastLayerEE)
328  delta_c = vecDeltas_[0];
329  else if (layer <= lastLayerFH)
330  delta_c = vecDeltas_[1];
331  else
332  delta_c = vecDeltas_[2];
333 
334  std::vector<size_t> rs =
335  sorted_indices(nd); // indices sorted by decreasing rho
336  std::vector<size_t> ds =
337  sort_by_delta(nd); // sort in decreasing distance to higher
338 
339  const unsigned int nd_size = nd.size();
340  for (unsigned int i = 0; i < nd_size; ++i) {
341 
342  if (nd[ds[i]].data.delta < delta_c)
343  break; // no more cluster centers to be looked at
344  if (dependSensor_) {
345 
346  float rho_c = kappa_ * nd[ds[i]].data.sigmaNoise;
347  if (nd[ds[i]].data.rho < rho_c)
348  continue; // set equal to kappa times noise threshold
349 
350  } else if (nd[ds[i]].data.rho * kappa_ < maxdensity)
351  continue;
352 
353  nd[ds[i]].data.clusterIndex = nClustersOnLayer;
354  if (verbosity_ < pINFO) {
355  std::cout << "Adding new cluster with index " << nClustersOnLayer
356  << std::endl;
357  std::cout << "Cluster center is hit " << ds[i] << std::endl;
358  }
359  nClustersOnLayer++;
360  }
361 
362  // at this point nClustersOnLayer is equal to the number of cluster centers -
363  // if it is zero we are done
364  if (nClustersOnLayer == 0)
365  return nClustersOnLayer;
366 
367  // assign remaining points to clusters, using the nearestHigher set from
368  // previous step (always set except
369  // for top density hit that is skipped...)
370  for (unsigned int oi = 1; oi < nd_size; ++oi) {
371  unsigned int i = rs[oi];
372  int ci = nd[i].data.clusterIndex;
373  if (ci ==
374  -1) { // clusterIndex is initialised with -1 if not yet used in cluster
375  nd[i].data.clusterIndex = nd[nd[i].data.nearestHigher].data.clusterIndex;
376  }
377  }
378 
379  // make room in the temporary cluster vector for the additional clusterIndex
380  // clusters
381  // from this layer
382  if (verbosity_ < pINFO) {
383  std::cout << "resizing cluster vector by " << nClustersOnLayer << std::endl;
384  }
385  clustersOnLayer.resize(nClustersOnLayer);
386 
387  // assign points closer than dc to other clusters to border region
388  // and find critical border density
389  std::vector<double> rho_b(nClustersOnLayer, 0.);
390  lp.clear();
391  lp.build(nd, bounds);
392  // now loop on all hits again :( and check: if there are hits from another
393  // cluster within d_c -> flag as border hit
394  for (unsigned int i = 0; i < nd_size; ++i) {
395  int ci = nd[i].data.clusterIndex;
396  bool flag_isolated = true;
397  if (ci != -1) {
398  KDTreeBox search_box(nd[i].dims[0] - delta_c, nd[i].dims[0] + delta_c,
399  nd[i].dims[1] - delta_c, nd[i].dims[1] + delta_c);
400  std::vector<KDNode> found;
401  lp.search(search_box, found);
402 
403  const unsigned int found_size = found.size();
404  for (unsigned int j = 0; j < found_size;
405  j++) { // start from 0 here instead of 1
406  // check if the hit is not within d_c of another cluster
407  if (found[j].data.clusterIndex != -1) {
408  float dist = distance(found[j].data, nd[i].data);
409  if (dist < delta_c && found[j].data.clusterIndex != ci) {
410  // in which case we assign it to the border
411  nd[i].data.isBorder = true;
412  break;
413  }
414  // because we are using two different containers, we have to make sure
415  // that we don't unflag the
416  // hit when it finds *itself* closer than delta_c
417  if (dist < delta_c && dist != 0. &&
418  found[j].data.clusterIndex == ci) {
419  // in this case it is not an isolated hit
420  // the dist!=0 is because the hit being looked at is also inside the
421  // search box and at dist==0
422  flag_isolated = false;
423  }
424  }
425  }
426  if (flag_isolated)
427  nd[i].data.isBorder =
428  true; // the hit is more than delta_c from any of its brethren
429  }
430  // check if this border hit has density larger than the current rho_b and
431  // update
432  if (nd[i].data.isBorder && rho_b[ci] < nd[i].data.rho)
433  rho_b[ci] = nd[i].data.rho;
434  } // end loop all hits
435 
436  // flag points in cluster with density < rho_b as halo points, then fill the
437  // cluster vector
438  for (unsigned int i = 0; i < nd_size; ++i) {
439  int ci = nd[i].data.clusterIndex;
440  if (ci != -1) {
441  if (nd[i].data.rho <= rho_b[ci])
442  nd[i].data.isHalo = true;
443  clustersOnLayer[ci].push_back(nd[i]);
444  if (verbosity_ < pINFO) {
445  std::cout << "Pushing hit " << i << " into cluster with index " << ci
446  << std::endl;
447  }
448  }
449  }
450 
451  // prepare the offset for the next layer if there is one
452  if (verbosity_ < pINFO) {
453  std::cout << "moving cluster offset by " << nClustersOnLayer << std::endl;
454  }
455  return nClustersOnLayer;
456 }
457 
458 // find local maxima within delta_c, marking the indices in the cluster
459 std::vector<unsigned>
460 HGCalImagingAlgo::findLocalMaximaInCluster(const std::vector<KDNode> &cluster) {
461  std::vector<unsigned> result;
462  std::vector<bool> seed(cluster.size(), true);
463  float delta_c = 2.;
464 
465  for (unsigned i = 0; i < cluster.size(); ++i) {
466  for (unsigned j = 0; j < cluster.size(); ++j) {
467  if (i != j and distance(cluster[i].data, cluster[j].data) < delta_c) {
468  if (cluster[i].data.weight < cluster[j].data.weight) {
469  seed[i] = false;
470  break;
471  }
472  }
473  }
474  }
475 
476  for (unsigned i = 0; i < cluster.size(); ++i) {
477  if (seed[i] && cluster[i].data.weight > 5e-4) {
478  // seed at i with energy cluster[i].weight
479  result.push_back(i);
480  }
481  }
482 
483  // Found result.size() sub-clusters in input cluster of length cluster.size()
484 
485  return result;
486 }
487 
489  const std::vector<KDNode> &hits, const std::vector<double> &fractions) {
490  double norm(0.0), x(0.0), y(0.0), z(0.0);
491  for (unsigned i = 0; i < hits.size(); ++i) {
492  const double weight = fractions[i] * hits[i].data.weight;
493  norm += weight;
494  x += weight * hits[i].data.x;
495  y += weight * hits[i].data.y;
496  z += weight * hits[i].data.z;
497  }
498  math::XYZPoint result(x, y, z);
499  result /= norm;
500  return result;
501 }
502 
504  const std::vector<KDNode> &hits, const std::vector<double> &fractions) {
505  double result = 0.0;
506  for (unsigned i = 0; i < hits.size(); ++i) {
507  result += fractions[i] * hits[i].data.weight;
508  }
509  return result;
510 }
511 
513  const std::vector<KDNode> &incluster, const std::vector<unsigned> &seeds,
514  std::vector<std::vector<double>> &outclusters) {
515  std::vector<bool> isaseed(incluster.size(), false);
516  outclusters.clear();
517  outclusters.resize(seeds.size());
518  std::vector<Point> centroids(seeds.size());
519  std::vector<double> energies(seeds.size());
520 
521  if (seeds.size() == 1) { // short circuit the case of a lone cluster
522  outclusters[0].clear();
523  outclusters[0].resize(incluster.size(), 1.0);
524  return;
525  }
526 
527  // saving seeds
528 
529  // create quick seed lookup
530  for (unsigned i = 0; i < seeds.size(); ++i) {
531  isaseed[seeds[i]] = true;
532  }
533 
534  // initialize clusters to be shared
535  // centroids start off at seed positions
536  // seeds always have fraction 1.0, to stabilize fit
537  // initializing fit
538  for (unsigned i = 0; i < seeds.size(); ++i) {
539  outclusters[i].resize(incluster.size(), 0.0);
540  for (unsigned j = 0; j < incluster.size(); ++j) {
541  if (j == seeds[i]) {
542  outclusters[i][j] = 1.0;
543  centroids[i] = math::XYZPoint(incluster[j].data.x, incluster[j].data.y,
544  incluster[j].data.z);
545  energies[i] = incluster[j].data.weight;
546  }
547  }
548  }
549 
550  // run the fit while we are less than max iterations, and clusters are still
551  // moving
552  const double minFracTot = 1e-20;
553  unsigned iter = 0;
554  const unsigned iterMax = 50;
556  const double stoppingTolerance = 1e-8;
557  const auto numberOfSeeds = seeds.size();
558  auto toleranceScaling =
559  numberOfSeeds > 2 ? (numberOfSeeds - 1) * (numberOfSeeds - 1) : 1;
560  std::vector<Point> prevCentroids;
561  std::vector<double> frac(numberOfSeeds), dist2(numberOfSeeds);
562  while (iter++ < iterMax && diff > stoppingTolerance * toleranceScaling) {
563  for (unsigned i = 0; i < incluster.size(); ++i) {
564  const Hexel &ihit = incluster[i].data;
565  double fracTot(0.0);
566  for (unsigned j = 0; j < numberOfSeeds; ++j) {
567  double fraction = 0.0;
568  double d2 = (std::pow(ihit.x - centroids[j].x(), 2.0) +
569  std::pow(ihit.y - centroids[j].y(), 2.0) +
570  std::pow(ihit.z - centroids[j].z(), 2.0)) /
571  sigma2_;
572  dist2[j] = d2;
573  // now we set the fractions up based on hit type
574  if (i == seeds[j]) { // this cluster's seed
575  fraction = 1.0;
576  } else if (isaseed[i]) {
577  fraction = 0.0;
578  } else {
579  fraction = energies[j] * std::exp(-0.5 * d2);
580  }
581  fracTot += fraction;
582  frac[j] = fraction;
583  }
584  // now that we have calculated all fractions for all hits
585  // assign the new fractions
586  for (unsigned j = 0; j < numberOfSeeds; ++j) {
587  if (fracTot > minFracTot || (i == seeds[j] && fracTot > 0.0)) {
588  outclusters[j][i] = frac[j] / fracTot;
589  } else {
590  outclusters[j][i] = 0.0;
591  }
592  }
593  }
594 
595  // save previous centroids
596  prevCentroids = std::move(centroids);
597  // finally update the position of the centroids from the last iteration
598  centroids.resize(numberOfSeeds);
599  double diff2 = 0.0;
600  for (unsigned i = 0; i < numberOfSeeds; ++i) {
601  centroids[i] = calculatePositionWithFraction(incluster, outclusters[i]);
602  energies[i] = calculateEnergyWithFraction(incluster, outclusters[i]);
603  // calculate convergence parameters
604  const double delta2 = (prevCentroids[i] - centroids[i]).perp2();
605  diff2 = std::max(delta2, diff2);
606  }
607  // update convergance parameter outside loop
608  diff = std::sqrt(diff2);
609  }
610 }
611 
613  // To support the TDR geometry and also the post-TDR one (v9 onwards), we
614  // need to change the logic of the vectors containing signal to noise and
615  // thresholds. The first 3 indices will keep on addressing the different
616  // thicknesses of the Silicon detectors, while the last one, number 3 (the
617  // fourth) will address the Scintillators. This change will support both
618  // geometries at the same time.
619 
620  if (initialized_)
621  return; // only need to calculate thresholds once
622 
623  initialized_ = true;
624 
625  std::vector<double> dummy;
626  const unsigned maxNumberOfThickIndices = 3;
627  dummy.resize(maxNumberOfThickIndices + 1, 0); // +1 to accomodate for the Scintillators
628  thresholds_.resize(maxlayer, dummy);
629  v_sigmaNoise_.resize(maxlayer, dummy);
630 
631  for (unsigned ilayer = 1; ilayer <= maxlayer; ++ilayer) {
632  for (unsigned ithick = 0; ithick < maxNumberOfThickIndices; ++ithick) {
633  float sigmaNoise =
634  0.001f * fcPerEle_ * nonAgedNoises_[ithick] * dEdXweights_[ilayer] /
635  (fcPerMip_[ithick] * thicknessCorrection_[ithick]);
636  thresholds_[ilayer - 1][ithick] = sigmaNoise * ecut_;
637  v_sigmaNoise_[ilayer - 1][ithick] = sigmaNoise;
638  }
639  float scintillators_sigmaNoise = 0.001f * noiseMip_ * dEdXweights_[ilayer];
640  thresholds_[ilayer - 1][maxNumberOfThickIndices] = ecut_ * scintillators_sigmaNoise;
641  v_sigmaNoise_[ilayer -1][maxNumberOfThickIndices] = scintillators_sigmaNoise;
642  }
643 
644 }
constexpr float energy() const
Definition: CaloRecHit.h:31
bool isHalfCell(const DetId &) const
Definition: RecHitTools.cc:348
double calculateEnergyWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
hgcal::RecHitTools rhtools_
math::XYZPoint calculatePositionWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
static const unsigned int maxlayer
double calculateLocalDensity(std::vector< KDNode > &, KDTree &, const unsigned int) const
int zside(const DetId &id) const
Definition: RecHitTools.cc:152
std::vector< std::vector< double > > thresholds_
void build(std::vector< KDTreeNodeInfoT< DATA, DIM > > &eltList, const KDTreeBoxT< DIM > &region)
constexpr const DetId & detid() const
Definition: CaloRecHit.h:35
std::vector< std::vector< double > > v_sigmaNoise_
math::XYZPoint Point
point in the space
Definition: weight.py:1
void search(const KDTreeBoxT< DIM > &searchBox, std::vector< KDTreeNodeInfoT< DATA, DIM > > &resRecHitList)
double distance2(const Hexel &pt1, const Hexel &pt2) const
void shareEnergy(const std::vector< KDNode > &, const std::vector< unsigned > &, std::vector< std::vector< double > > &)
static const unsigned int lastLayerEE
double distance(const Hexel &pt1, const Hexel &pt2) const
T sqrt(T t)
Definition: SSEVec.h:18
static const unsigned int lastLayerFH
std::vector< reco::BasicCluster > clusters_v_
reco::CaloCluster::AlgoId algoId_
std::vector< std::vector< KDNode > > points_
std::vector< double > vecDeltas_
std::vector< double > dEdXweights_
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< size_t > sorted_indices(const std::vector< T > &v)
std::vector< std::vector< std::vector< KDNode > > > layerClustersPerLayer_
std::vector< reco::BasicCluster > getClusters(bool)
std::vector< unsigned > findLocalMaximaInCluster(const std::vector< KDNode > &)
std::vector< std::array< float, 2 > > minpos_
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:166
Definition: DetId.h:18
std::vector< double > thicknessCorrection_
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:304
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:189
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
T perp2() const
Squared magnitude of transverse component.
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:133
std::vector< double > fcPerMip_
std::vector< double > nonAgedNoises_
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
size_type size() const
double calculateDistanceToHigher(std::vector< KDNode > &) const
static int position[264][3]
Definition: ReadPGInfo.cc:509
VerbosityLevel verbosity_
math::XYZPoint calculatePosition(std::vector< KDNode > &) const
std::vector< std::array< float, 2 > > maxpos_
void populate(const HGCRecHitCollection &hits)
std::vector< size_t > sort_by_delta(const std::vector< KDNode > &v) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:511
int findAndAssignClusters(std::vector< KDNode > &, KDTree &, double, KDTreeBox &, const unsigned int, std::vector< std::vector< KDNode > > &) const