CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes
HGCalImagingAlgo Class Reference

#include <HGCalImagingAlgo.h>

Inheritance diagram for HGCalImagingAlgo:
HGCalClusteringAlgoBase

Classes

struct  Hexel
 

Public Types

typedef math::XYZPoint Point
 point in the space More...
 
- Public Types inherited from HGCalClusteringAlgoBase
enum  VerbosityLevel { pDEBUG = 0, pWARNING = 1, pINFO = 2, pERROR = 3 }
 

Public Member Functions

void computeThreshold ()
 
std::vector< reco::BasicClustergetClusters (bool) override
 
Density getDensity () override
 
void getEventSetupPerAlgorithm (const edm::EventSetup &es) override
 
 HGCalImagingAlgo (const edm::ParameterSet &ps, edm::ConsumesCollector iC)
 
void makeClusters () override
 
void populate (const HGCRecHitCollection &hits) override
 
void reset () override
 
 ~HGCalImagingAlgo () override
 
- Public Member Functions inherited from HGCalClusteringAlgoBase
void getEventSetup (const edm::EventSetup &es)
 
 HGCalClusteringAlgoBase (VerbosityLevel v, reco::CaloCluster::AlgoId algo, edm::ConsumesCollector iC)
 
void setAlgoId (reco::CaloCluster::AlgoId algo, bool isNose=false)
 
void setVerbosity (VerbosityLevel the_verbosity)
 
virtual ~HGCalClusteringAlgoBase ()
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &iDesc)
 

Private Types

typedef KDTreeNodeInfo< HexelKDNode
 
typedef KDTreeLinkerAlgo< HexelKDTree
 

Private Member Functions

double calculateDistanceToHigher (std::vector< KDNode > &) const
 
double calculateEnergyWithFraction (const std::vector< KDNode > &, const std::vector< double > &)
 
double calculateLocalDensity (std::vector< KDNode > &, KDTree &, const unsigned int) const
 
math::XYZPoint calculatePosition (std::vector< KDNode > &) const
 
math::XYZPoint calculatePositionWithFraction (const std::vector< KDNode > &, const std::vector< double > &)
 
double distance (const Hexel &pt1, const Hexel &pt2) const
 
double distance2 (const Hexel &pt1, const Hexel &pt2) const
 
int findAndAssignClusters (std::vector< KDNode > &, KDTree &, double, KDTreeBox< 2 > &, const unsigned int, std::vector< std::vector< KDNode >> &) const
 
std::vector< unsigned > findLocalMaximaInCluster (const std::vector< KDNode > &)
 
void setDensity (const std::vector< KDNode > &nd)
 
void shareEnergy (const std::vector< KDNode > &, const std::vector< unsigned > &, std::vector< std::vector< double >> &)
 
std::vector< size_t > sort_by_delta (const std::vector< KDNode > &v) const
 

Private Attributes

std::vector< reco::BasicClusterclusters_v_
 
std::vector< double > dEdXweights_
 
Density density_
 
bool dependSensor_
 
double ecut_
 
double fcPerEle_
 
std::vector< double > fcPerMip_
 
bool initialized_
 
double kappa_
 
std::vector< std::vector< std::vector< KDNode > > > layerClustersPerLayer_
 
std::vector< std::array< float, 2 > > maxpos_
 
std::vector< std::array< float, 2 > > minpos_
 
double noiseMip_
 
std::vector< double > nonAgedNoises_
 
std::vector< std::vector< KDNode > > points_
 
std::vector< double > positionDeltaRho_c_
 
double sigma2_
 
std::vector< std::vector< double > > sigmaNoise_
 
std::vector< double > thicknessCorrection_
 
std::vector< std::vector< double > > thresholds_
 
std::vector< double > thresholdW0_
 
std::vector< double > vecDeltas_
 

Additional Inherited Members

- Public Attributes inherited from HGCalClusteringAlgoBase
unsigned int firstLayerBH_
 
bool isNose_
 
unsigned int lastLayerEE_
 
unsigned int lastLayerFH_
 
unsigned int maxlayer_
 
int scintMaxIphi_
 
- Protected Attributes inherited from HGCalClusteringAlgoBase
reco::CaloCluster::AlgoId algoId_
 
edm::ESGetToken< CaloGeometry, CaloGeometryRecordcaloGeomToken_
 
std::vector< reco::BasicClusterclusters_v_
 
hgcal::RecHitTools rhtools_
 
VerbosityLevel verbosity_
 

Detailed Description

Definition at line 31 of file HGCalImagingAlgo.h.

Member Typedef Documentation

◆ KDNode

Definition at line 221 of file HGCalImagingAlgo.h.

◆ KDTree

Definition at line 220 of file HGCalImagingAlgo.h.

◆ Point

point in the space

Definition at line 120 of file HGCalImagingAlgo.h.

Constructor & Destructor Documentation

◆ HGCalImagingAlgo()

HGCalImagingAlgo::HGCalImagingAlgo ( const edm::ParameterSet ps,
edm::ConsumesCollector  iC 
)
inline

Definition at line 33 of file HGCalImagingAlgo.h.

35  (HGCalClusteringAlgoBase::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity", 3),
37  iC),
38  thresholdW0_(ps.getParameter<std::vector<double>>("thresholdW0")),
39  positionDeltaRho_c_(ps.getParameter<std::vector<double>>("positionDeltaRho_c")),
40  vecDeltas_(ps.getParameter<std::vector<double>>("deltac")),
41  kappa_(ps.getParameter<double>("kappa")),
42  ecut_(ps.getParameter<double>("ecut")),
43  sigma2_(1.0),
44  dependSensor_(ps.getParameter<bool>("dependSensor")),
45  dEdXweights_(ps.getParameter<std::vector<double>>("dEdXweights")),
46  thicknessCorrection_(ps.getParameter<std::vector<double>>("thicknessCorrection")),
47  fcPerMip_(ps.getParameter<std::vector<double>>("fcPerMip")),
48  fcPerEle_(ps.getParameter<double>("fcPerEle")),
49  nonAgedNoises_(ps.getParameter<edm::ParameterSet>("noises").getParameter<std::vector<double>>("values")),
50  noiseMip_(ps.getParameter<edm::ParameterSet>("noiseMip").getParameter<double>("noise_MIP")),
51  initialized_(false) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HGCalClusteringAlgoBase(VerbosityLevel v, reco::CaloCluster::AlgoId algo, edm::ConsumesCollector iC)
std::vector< double > positionDeltaRho_c_
std::vector< double > thresholdW0_
T getUntrackedParameter(std::string const &, T const &) const
std::vector< double > vecDeltas_
std::vector< double > dEdXweights_
std::vector< double > thicknessCorrection_
std::vector< double > fcPerMip_
std::vector< double > nonAgedNoises_

◆ ~HGCalImagingAlgo()

HGCalImagingAlgo::~HGCalImagingAlgo ( )
inlineoverride

Definition at line 53 of file HGCalImagingAlgo.h.

53 {}

Member Function Documentation

◆ calculateDistanceToHigher()

double HGCalImagingAlgo::calculateDistanceToHigher ( std::vector< KDNode > &  nd) const
private

Definition at line 291 of file HGCalImagingAlgo.cc.

References data, cms::cuda::for(), mps_fire::i, dqmiolumiharvest::j, hgcal_clustering::sorted_indices(), mathSSE::sqrt(), and createJobs::tmp.

291  {
292  // sort vector of Hexels by decreasing local density
293  std::vector<size_t> rs = sorted_indices(nd);
294 
295  double maxdensity = 0.0;
296  int nearestHigher = -1;
297 
298  if (!rs.empty())
299  maxdensity = nd[rs[0]].data.rho;
300  else
301  return maxdensity; // there are no hits
302  double dist2 = 0.;
303  // start by setting delta for the highest density hit to
304  // the most distant hit - this is a convention
305 
306  for (auto &j : nd) {
307  double tmp = distance2(nd[rs[0]].data, j.data);
308  if (tmp > dist2)
309  dist2 = tmp;
310  }
311  nd[rs[0]].data.delta = std::sqrt(dist2);
312  nd[rs[0]].data.nearestHigher = nearestHigher;
313 
314  // now we save the largest distance as a starting point
315  const double max_dist2 = dist2;
316  const unsigned int nd_size = nd.size();
317 
318  for (unsigned int oi = 1; oi < nd_size; ++oi) { // start from second-highest density
319  dist2 = max_dist2;
320  unsigned int i = rs[oi];
321  // we only need to check up to oi since hits
322  // are ordered by decreasing density
323  // and all points coming BEFORE oi are guaranteed to have higher rho
324  // and the ones AFTER to have lower rho
325  for (unsigned int oj = 0; oj < oi; ++oj) {
326  unsigned int j = rs[oj];
327  // Limit the search box
328  if ((nd[i].data.x - nd[j].data.x) * (nd[i].data.x - nd[j].data.x) > dist2)
329  continue;
330  if ((nd[i].data.y - nd[j].data.y) * (nd[i].data.y - nd[j].data.y) > dist2)
331  continue;
332  double tmp = distance2(nd[i].data, nd[j].data);
333  if (tmp <= dist2) { // this "<=" instead of "<" addresses the (rare) case
334  // when there are only two hits
335  dist2 = tmp;
336  nearestHigher = j;
337  }
338  }
339  nd[i].data.delta = std::sqrt(dist2);
340  nd[i].data.nearestHigher = nearestHigher; // this uses the original unsorted hitlist
341  }
342  return maxdensity;
343 }
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
double distance2(const Hexel &pt1, const Hexel &pt2) const
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< size_t > sorted_indices(const std::vector< T > &v)
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
tmp
align.sh
Definition: createJobs.py:716

◆ calculateEnergyWithFraction()

double HGCalImagingAlgo::calculateEnergyWithFraction ( const std::vector< KDNode > &  hits,
const std::vector< double > &  fractions 
)
private

Definition at line 522 of file HGCalImagingAlgo.cc.

References hfClusterShapes_cfi::hits, mps_fire::i, and mps_fire::result.

523  {
524  double result = 0.0;
525  for (unsigned i = 0; i < hits.size(); ++i) {
526  result += fractions[i] * hits[i].data.weight;
527  }
528  return result;
529 }

◆ calculateLocalDensity()

double HGCalImagingAlgo::calculateLocalDensity ( std::vector< KDNode > &  nd,
KDTree lp,
const unsigned int  layer 
) const
private

Definition at line 262 of file HGCalImagingAlgo.cc.

References data, HLT_2022v15_cff::distance, newFWLiteAna::found, mps_fire::i, dqmiolumiharvest::j, SiStripPI::max, and KDTreeLinkerAlgo< DATA, DIM >::search().

262  {
263  double maxdensity = 0.;
264  float delta_c; // maximum search distance (critical distance) for local
265  // density calculation
266  if (layer <= lastLayerEE_)
267  delta_c = vecDeltas_[0];
268  else if (layer < firstLayerBH_)
269  delta_c = vecDeltas_[1];
270  else
271  delta_c = vecDeltas_[2];
272 
273  // for each node calculate local density rho and store it
274  for (unsigned int i = 0; i < nd.size(); ++i) {
275  // speec up search by looking within +/- delta_c window only
276  KDTreeBox search_box(
277  nd[i].dims[0] - delta_c, nd[i].dims[0] + delta_c, nd[i].dims[1] - delta_c, nd[i].dims[1] + delta_c);
278  std::vector<Hexel> found;
279  lp.search(search_box, found);
280  const unsigned int found_size = found.size();
281  for (unsigned int j = 0; j < found_size; j++) {
282  if (distance(nd[i].data, found[j]) < delta_c) {
283  nd[i].data.rho += found[j].weight;
284  maxdensity = std::max(maxdensity, nd[i].data.rho);
285  }
286  } // end loop found
287  } // end loop nodes
288  return maxdensity;
289 }
double distance(const Hexel &pt1, const Hexel &pt2) const
constexpr std::array< uint8_t, layerIndexSize > layer
std::vector< double > vecDeltas_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79

◆ calculatePosition()

math::XYZPoint HGCalImagingAlgo::calculatePosition ( std::vector< KDNode > &  v) const
private

Definition at line 196 of file HGCalImagingAlgo.cc.

References data, mps_fire::i, heavyIonCSV_trainingSettings::idx, dqm-mbProfile::log, SiStripPI::max, findQualityFiles::v, and x.

196  {
197  float total_weight = 0.f;
198  float x = 0.f;
199  float y = 0.f;
200 
201  unsigned int v_size = v.size();
202  unsigned int maxEnergyIndex = 0;
203  float maxEnergyValue = 0;
204 
205  // loop over hits in cluster candidate
206  // determining the maximum energy hit
207  for (unsigned int i = 0; i < v_size; i++) {
208  if (v[i].data.weight > maxEnergyValue) {
209  maxEnergyValue = v[i].data.weight;
210  maxEnergyIndex = i;
211  }
212  }
213 
214  // Si cell or Scintillator. Used to set approach and parameters
215  int thick = rhtools_.getSiThickIndex(v[maxEnergyIndex].data.detid);
216 
217  // for hits within positionDeltaRho_c_ from maximum energy hit
218  // build up weight for energy-weighted position
219  // and save corresponding hits indices
220  std::vector<unsigned int> innerIndices;
221  for (unsigned int i = 0; i < v_size; i++) {
222  if (thick == -1 || distance2(v[i].data, v[maxEnergyIndex].data) < positionDeltaRho_c_[thick]) {
223  innerIndices.push_back(i);
224 
225  float rhEnergy = v[i].data.weight;
226  total_weight += rhEnergy;
227  // just fill x, y for scintillator
228  // for Si it is overwritten later anyway
229  if (thick == -1) {
230  x += v[i].data.x * rhEnergy;
231  y += v[i].data.y * rhEnergy;
232  }
233  }
234  }
235  // just loop on reduced vector of interesting indices
236  // to compute log weighting
237  if (thick != -1 && total_weight != 0.) { // Silicon case
238  float total_weight_log = 0.f;
239  float x_log = 0.f;
240  float y_log = 0.f;
241  for (auto idx : innerIndices) {
242  float rhEnergy = v[idx].data.weight;
243  if (rhEnergy == 0.)
244  continue;
245  float Wi = std::max(thresholdW0_[thick] + log(rhEnergy / total_weight), 0.);
246  x_log += v[idx].data.x * Wi;
247  y_log += v[idx].data.y * Wi;
248  total_weight_log += Wi;
249  }
250  total_weight = total_weight_log;
251  x = x_log;
252  y = y_log;
253  }
254 
255  if (total_weight != 0.) {
256  auto inv_tot_weight = 1. / total_weight;
257  return math::XYZPoint(x * inv_tot_weight, y * inv_tot_weight, v[maxEnergyIndex].data.z);
258  }
259  return math::XYZPoint(0, 0, 0);
260 }
double distance2(const Hexel &pt1, const Hexel &pt2) const
std::vector< double > positionDeltaRho_c_
std::vector< double > thresholdW0_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:205

◆ calculatePositionWithFraction()

math::XYZPoint HGCalImagingAlgo::calculatePositionWithFraction ( const std::vector< KDNode > &  hits,
const std::vector< double > &  fractions 
)
private

Definition at line 507 of file HGCalImagingAlgo.cc.

References hfClusterShapes_cfi::hits, mps_fire::i, mps_fire::result, mps_merge::weight, and x.

508  {
509  double norm(0.0), x(0.0), y(0.0), z(0.0);
510  for (unsigned i = 0; i < hits.size(); ++i) {
511  const double weight = fractions[i] * hits[i].data.weight;
512  norm += weight;
513  x += weight * hits[i].data.x;
514  y += weight * hits[i].data.y;
515  z += weight * hits[i].data.z;
516  }
518  result /= norm;
519  return result;
520 }
Definition: weight.py:1
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12

◆ computeThreshold()

void HGCalImagingAlgo::computeThreshold ( )

Definition at line 628 of file HGCalImagingAlgo.cc.

References hgcalLayerClusters_cfi::maxNumberOfThickIndices.

628  {
629  // To support the TDR geometry and also the post-TDR one (v9 onwards), we
630  // need to change the logic of the vectors containing signal to noise and
631  // thresholds. The first 3 indices will keep on addressing the different
632  // thicknesses of the Silicon detectors, while the last one, number 3 (the
633  // fourth) will address the Scintillators. This change will support both
634  // geometries at the same time.
635 
636  if (initialized_)
637  return; // only need to calculate thresholds once
638 
639  initialized_ = true;
640 
641  std::vector<double> dummy;
642  const unsigned maxNumberOfThickIndices = 3;
643  dummy.resize(maxNumberOfThickIndices + 1, 0); // +1 to accomodate for the Scintillators
644  thresholds_.resize(maxlayer_, dummy);
645  sigmaNoise_.resize(maxlayer_, dummy);
646 
647  for (unsigned ilayer = 1; ilayer <= maxlayer_; ++ilayer) {
648  for (unsigned ithick = 0; ithick < maxNumberOfThickIndices; ++ithick) {
649  float sigmaNoise = 0.001f * fcPerEle_ * nonAgedNoises_[ithick] * dEdXweights_[ilayer] /
650  (fcPerMip_[ithick] * thicknessCorrection_[ithick]);
651  thresholds_[ilayer - 1][ithick] = sigmaNoise * ecut_;
652  sigmaNoise_[ilayer - 1][ithick] = sigmaNoise;
653  }
654  float scintillators_sigmaNoise = 0.001f * noiseMip_ * dEdXweights_[ilayer];
655  thresholds_[ilayer - 1][maxNumberOfThickIndices] = ecut_ * scintillators_sigmaNoise;
656  sigmaNoise_[ilayer - 1][maxNumberOfThickIndices] = scintillators_sigmaNoise;
657  }
658 }
std::vector< std::vector< double > > thresholds_
std::vector< double > dEdXweights_
std::vector< std::vector< double > > sigmaNoise_
std::vector< double > thicknessCorrection_
std::vector< double > fcPerMip_
std::vector< double > nonAgedNoises_

◆ distance()

double HGCalImagingAlgo::distance ( const Hexel pt1,
const Hexel pt2 
) const
inlineprivate

Definition at line 244 of file HGCalImagingAlgo.h.

References distance2(), HLT_2022v15_cff::pt1, HLT_2022v15_cff::pt2, and mathSSE::sqrt().

244  { //2-d distance on the layer (x-y)
245  return std::sqrt(distance2(pt1, pt2));
246  }
double distance2(const Hexel &pt1, const Hexel &pt2) const
T sqrt(T t)
Definition: SSEVec.h:19

◆ distance2()

double HGCalImagingAlgo::distance2 ( const Hexel pt1,
const Hexel pt2 
) const
inlineprivate

Definition at line 239 of file HGCalImagingAlgo.h.

References PVValHelper::dx, PVValHelper::dy, HLT_2022v15_cff::pt1, and HLT_2022v15_cff::pt2.

Referenced by distance().

239  { //distance squared
240  const double dx = pt1.x - pt2.x;
241  const double dy = pt1.y - pt2.y;
242  return (dx * dx + dy * dy);
243  } //distance squaredq

◆ fillPSetDescription()

static void HGCalImagingAlgo::fillPSetDescription ( edm::ParameterSetDescription iDesc)
inlinestatic

Definition at line 89 of file HGCalImagingAlgo.h.

References edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addUntracked(), and AlCaHLTBitMon_QueryRunRegistry::string.

89  {
90  iDesc.add<std::vector<double>>("thresholdW0", {2.9, 2.9, 2.9});
91  iDesc.add<std::vector<double>>("positionDeltaRho_c", {1.3, 1.3, 1.3});
92  iDesc.add<std::vector<double>>("deltac",
93  {
94  2.0,
95  2.0,
96  5.0,
97  });
98  iDesc.add<bool>("dependSensor", true);
99  iDesc.add<double>("ecut", 3.0);
100  iDesc.add<double>("kappa", 9.0);
101  iDesc.addUntracked<unsigned int>("verbosity", 3);
102  iDesc.add<std::vector<double>>("dEdXweights", {});
103  iDesc.add<std::vector<double>>("thicknessCorrection", {});
104  iDesc.add<std::vector<double>>("fcPerMip", {});
105  iDesc.add<double>("fcPerEle", 0.0);
106  edm::ParameterSetDescription descNestedNoises;
107  descNestedNoises.add<std::vector<double>>("values", {});
108  iDesc.add<edm::ParameterSetDescription>("noises", descNestedNoises);
109  edm::ParameterSetDescription descNestedNoiseMIP;
110  descNestedNoiseMIP.add<bool>("scaleByDose", false);
111  descNestedNoiseMIP.add<double>("scaleByDoseFactor", 1.);
112  iDesc.add<edm::ParameterSetDescription>("scaleByDose", descNestedNoiseMIP);
113  descNestedNoiseMIP.add<std::string>("doseMap", "");
114  iDesc.add<edm::ParameterSetDescription>("doseMap", descNestedNoiseMIP);
115  descNestedNoiseMIP.add<double>("noise_MIP", 1. / 100.);
116  iDesc.add<edm::ParameterSetDescription>("noiseMip", descNestedNoiseMIP);
117  }
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ findAndAssignClusters()

int HGCalImagingAlgo::findAndAssignClusters ( std::vector< KDNode > &  nd,
KDTree lp,
double  maxdensity,
KDTreeBox< 2 > &  bounds,
const unsigned int  layer,
std::vector< std::vector< KDNode >> &  clustersOnLayer 
) const
private

Definition at line 344 of file HGCalImagingAlgo.cc.

References KDTreeLinkerAlgo< DATA, DIM >::build(), KDTreeLinkerAlgo< DATA, DIM >::clear(), gather_cfg::cout, data, HLT_2022v15_cff::distance, newFWLiteAna::found, mps_fire::i, dqmiolumiharvest::j, KDTreeLinkerAlgo< DATA, DIM >::search(), and hgcal_clustering::sorted_indices().

349  {
350  // this is called once per layer and endcap...
351  // so when filling the cluster temporary vector of Hexels we resize each time
352  // by the number of clusters found. This is always equal to the number of
353  // cluster centers...
354 
355  unsigned int nClustersOnLayer = 0;
356  float delta_c; // critical distance
357  if (layer <= lastLayerEE_)
358  delta_c = vecDeltas_[0];
359  else if (layer < firstLayerBH_)
360  delta_c = vecDeltas_[1];
361  else
362  delta_c = vecDeltas_[2];
363 
364  std::vector<size_t> rs = sorted_indices(nd); // indices sorted by decreasing rho
365  std::vector<size_t> ds = sort_by_delta(nd); // sort in decreasing distance to higher
366 
367  const unsigned int nd_size = nd.size();
368  for (unsigned int i = 0; i < nd_size; ++i) {
369  if (nd[ds[i]].data.delta < delta_c)
370  break; // no more cluster centers to be looked at
371  if (dependSensor_) {
372  float rho_c = kappa_ * nd[ds[i]].data.sigmaNoise;
373  if (nd[ds[i]].data.rho < rho_c)
374  continue; // set equal to kappa times noise threshold
375 
376  } else if (nd[ds[i]].data.rho * kappa_ < maxdensity)
377  continue;
378 
379  nd[ds[i]].data.clusterIndex = nClustersOnLayer;
380  if (verbosity_ < pINFO) {
381  std::cout << "Adding new cluster with index " << nClustersOnLayer << std::endl;
382  std::cout << "Cluster center is hit " << ds[i] << std::endl;
383  }
384  nClustersOnLayer++;
385  }
386 
387  // at this point nClustersOnLayer is equal to the number of cluster centers -
388  // if it is zero we are done
389  if (nClustersOnLayer == 0)
390  return nClustersOnLayer;
391 
392  // assign remaining points to clusters, using the nearestHigher set from
393  // previous step (always set except
394  // for top density hit that is skipped...)
395  for (unsigned int oi = 1; oi < nd_size; ++oi) {
396  unsigned int i = rs[oi];
397  int ci = nd[i].data.clusterIndex;
398  if (ci == -1) { // clusterIndex is initialised with -1 if not yet used in cluster
399  nd[i].data.clusterIndex = nd[nd[i].data.nearestHigher].data.clusterIndex;
400  }
401  }
402 
403  // make room in the temporary cluster vector for the additional clusterIndex
404  // clusters
405  // from this layer
406  if (verbosity_ < pINFO) {
407  std::cout << "resizing cluster vector by " << nClustersOnLayer << std::endl;
408  }
409  clustersOnLayer.resize(nClustersOnLayer);
410 
411  // assign points closer than dc to other clusters to border region
412  // and find critical border density
413  std::vector<double> rho_b(nClustersOnLayer, 0.);
414  lp.clear();
415  lp.build(nd, bounds);
416  // now loop on all hits again :( and check: if there are hits from another
417  // cluster within d_c -> flag as border hit
418  for (unsigned int i = 0; i < nd_size; ++i) {
419  int ci = nd[i].data.clusterIndex;
420  bool flag_isolated = true;
421  if (ci != -1) {
422  KDTreeBox search_box(
423  nd[i].dims[0] - delta_c, nd[i].dims[0] + delta_c, nd[i].dims[1] - delta_c, nd[i].dims[1] + delta_c);
424  std::vector<Hexel> found;
425  lp.search(search_box, found);
426 
427  const unsigned int found_size = found.size();
428  for (unsigned int j = 0; j < found_size; j++) { // start from 0 here instead of 1
429  // check if the hit is not within d_c of another cluster
430  if (found[j].clusterIndex != -1) {
431  float dist = distance(found[j], nd[i].data);
432  if (dist < delta_c && found[j].clusterIndex != ci) {
433  // in which case we assign it to the border
434  nd[i].data.isBorder = true;
435  break;
436  }
437  // because we are using two different containers, we have to make sure
438  // that we don't unflag the
439  // hit when it finds *itself* closer than delta_c
440  if (dist < delta_c && dist != 0. && found[j].clusterIndex == ci) {
441  // in this case it is not an isolated hit
442  // the dist!=0 is because the hit being looked at is also inside the
443  // search box and at dist==0
444  flag_isolated = false;
445  }
446  }
447  }
448  if (flag_isolated)
449  nd[i].data.isBorder = true; // the hit is more than delta_c from any of its brethren
450  }
451  // check if this border hit has density larger than the current rho_b and
452  // update
453  if (nd[i].data.isBorder && rho_b[ci] < nd[i].data.rho)
454  rho_b[ci] = nd[i].data.rho;
455  } // end loop all hits
456 
457  // flag points in cluster with density < rho_b as halo points, then fill the
458  // cluster vector
459  for (unsigned int i = 0; i < nd_size; ++i) {
460  int ci = nd[i].data.clusterIndex;
461  if (ci != -1) {
462  if (nd[i].data.rho <= rho_b[ci])
463  nd[i].data.isHalo = true;
464  clustersOnLayer[ci].push_back(nd[i]);
465  if (verbosity_ < pINFO) {
466  std::cout << "Pushing hit " << i << " into cluster with index " << ci << std::endl;
467  }
468  }
469  }
470 
471  // prepare the offset for the next layer if there is one
472  if (verbosity_ < pINFO) {
473  std::cout << "moving cluster offset by " << nClustersOnLayer << std::endl;
474  }
475  return nClustersOnLayer;
476 }
double distance(const Hexel &pt1, const Hexel &pt2) const
constexpr std::array< uint8_t, layerIndexSize > layer
std::vector< size_t > sort_by_delta(const std::vector< KDNode > &v) const
std::vector< size_t > sorted_indices(const std::vector< T > &v)
std::vector< double > vecDeltas_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79

◆ findLocalMaximaInCluster()

std::vector< unsigned > HGCalImagingAlgo::findLocalMaximaInCluster ( const std::vector< KDNode > &  cluster)
private

Definition at line 479 of file HGCalImagingAlgo.cc.

References data, HLT_2022v15_cff::distance, MillePedeFileConverter_cfg::e, mps_fire::i, dqmiolumiharvest::j, mps_fire::result, and fileCollector::seed.

479  {
480  std::vector<unsigned> result;
481  std::vector<bool> seed(cluster.size(), true);
482  float delta_c = 2.;
483 
484  for (unsigned i = 0; i < cluster.size(); ++i) {
485  for (unsigned j = 0; j < cluster.size(); ++j) {
486  if (i != j and distance(cluster[i].data, cluster[j].data) < delta_c) {
487  if (cluster[i].data.weight < cluster[j].data.weight) {
488  seed[i] = false;
489  break;
490  }
491  }
492  }
493  }
494 
495  for (unsigned i = 0; i < cluster.size(); ++i) {
496  if (seed[i] && cluster[i].data.weight > 5e-4) {
497  // seed at i with energy cluster[i].weight
498  result.push_back(i);
499  }
500  }
501 
502  // Found result.size() sub-clusters in input cluster of length cluster.size()
503 
504  return result;
505 }
double distance(const Hexel &pt1, const Hexel &pt2) const
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79

◆ getClusters()

std::vector< reco::BasicCluster > HGCalImagingAlgo::getClusters ( bool  doSharing)
overridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 119 of file HGCalImagingAlgo.cc.

References AlignmentPI::calculatePosition(), gather_cfg::cout, data, reco::CaloID::DET_HGCAL_ENDCAP, hgcalLayerClusters_cfi::doSharing, MillePedeFileConverter_cfg::e, HCALHighEnergyHPDFilter_cfi::energy, HLT_2022v15_cff::fraction, mps_fire::i, hgcal_clustering::max_index(), position, and DetachedQuadStep_cff::seeds.

119  {
121  std::vector<std::pair<DetId, float>> thisCluster;
122  for (auto &clsOnLayer : layerClustersPerLayer_) {
123  for (unsigned int i = 0; i < clsOnLayer.size(); ++i) {
124  double energy = 0;
125  Point position;
126 
127  //Will save the maximum density hit of the cluster
128  size_t rsmax = max_index(clsOnLayer[i]);
129 
130  if (doSharing) {
131  std::vector<unsigned> seeds = findLocalMaximaInCluster(clsOnLayer[i]);
132  // sharing found seeds.size() sub-cluster seeds in cluster i
133 
134  std::vector<std::vector<double>> fractions;
135  // first pass can have noise it in
136  shareEnergy(clsOnLayer[i], seeds, fractions);
137 
138  // reset and run second pass after vetoing seeds
139  // that result in trivial clusters (less than 2 effective cells)
140 
141  for (unsigned isub = 0; isub < fractions.size(); ++isub) {
142  double effective_hits = 0.0;
143  double energy = calculateEnergyWithFraction(clsOnLayer[i], fractions[isub]);
144  Point position = calculatePositionWithFraction(clsOnLayer[i], fractions[isub]);
145 
146  for (unsigned ihit = 0; ihit < fractions[isub].size(); ++ihit) {
147  const double fraction = fractions[isub][ihit];
148  if (fraction > 1e-7) {
149  effective_hits += fraction;
150  thisCluster.emplace_back(clsOnLayer[i][ihit].data.detid, fraction);
151  }
152  }
153 
154  if (verbosity_ < pINFO) {
155  std::cout << "\t******** NEW CLUSTER (SHARING) ********" << std::endl;
156  std::cout << "\tEff. No. of cells = " << effective_hits << std::endl;
157  std::cout << "\t Energy = " << energy << std::endl;
158  std::cout << "\t Phi = " << position.phi() << std::endl;
159  std::cout << "\t Eta = " << position.eta() << std::endl;
160  std::cout << "\t*****************************" << std::endl;
161  }
162  clusters_v_.emplace_back(energy, position, caloID, thisCluster, algoId_);
163  if (!clusters_v_.empty()) {
164  clusters_v_.back().setSeed(clsOnLayer[i][rsmax].data.detid);
165  }
166  thisCluster.clear();
167  }
168  } else {
169  position = calculatePosition(clsOnLayer[i]); // energy-weighted position
170  // std::vector< KDNode >::iterator it;
171  for (auto &it : clsOnLayer[i]) {
172  energy += it.data.isHalo ? 0. : it.data.weight;
173  // use fraction to store whether this is a Halo hit or not
174  thisCluster.emplace_back(it.data.detid, (it.data.isHalo ? 0.f : 1.f));
175  }
176  if (verbosity_ < pINFO) {
177  std::cout << "******** NEW CLUSTER (HGCIA) ********" << std::endl;
178  std::cout << "Index " << i << std::endl;
179  std::cout << "No. of cells = " << clsOnLayer[i].size() << std::endl;
180  std::cout << " Energy = " << energy << std::endl;
181  std::cout << " Phi = " << position.phi() << std::endl;
182  std::cout << " Eta = " << position.eta() << std::endl;
183  std::cout << "*****************************" << std::endl;
184  }
185  clusters_v_.emplace_back(energy, position, caloID, thisCluster, algoId_);
186  if (!clusters_v_.empty()) {
187  clusters_v_.back().setSeed(clsOnLayer[i][rsmax].data.detid);
188  }
189  thisCluster.clear();
190  }
191  }
192  }
193  return clusters_v_;
194 }
math::XYZPoint calculatePosition(std::vector< KDNode > &) const
double calculateEnergyWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
math::XYZPoint calculatePositionWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
void shareEnergy(const std::vector< KDNode > &, const std::vector< unsigned > &, std::vector< std::vector< double >> &)
std::vector< reco::BasicCluster > clusters_v_
reco::CaloCluster::AlgoId algoId_
std::vector< unsigned > findLocalMaximaInCluster(const std::vector< KDNode > &)
std::vector< std::vector< std::vector< KDNode > > > layerClustersPerLayer_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
Structure Point Contains parameters of Gaussian fits to DMRs.
static int position[264][3]
Definition: ReadPGInfo.cc:289
size_t max_index(const std::vector< T > &v)

◆ getDensity()

Density HGCalImagingAlgo::getDensity ( )
overridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 668 of file HGCalImagingAlgo.cc.

668 { return density_; }

◆ getEventSetupPerAlgorithm()

void HGCalImagingAlgo::getEventSetupPerAlgorithm ( const edm::EventSetup es)
overridevirtual

Reimplemented from HGCalClusteringAlgoBase.

Definition at line 17 of file HGCalImagingAlgo.cc.

17  {
18  points_.clear();
19  minpos_.clear();
20  maxpos_.clear();
21  points_.resize(2 * (maxlayer_ + 1));
22  minpos_.resize(2 * (maxlayer_ + 1), {{0.0f, 0.0f}});
23  maxpos_.resize(2 * (maxlayer_ + 1), {{0.0f, 0.0f}});
24 }
std::vector< std::array< float, 2 > > minpos_
std::vector< std::array< float, 2 > > maxpos_
std::vector< std::vector< KDNode > > points_

◆ makeClusters()

void HGCalImagingAlgo::makeClusters ( )
overridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 93 of file HGCalImagingAlgo.cc.

References KDTreeLinkerAlgo< DATA, DIM >::build(), mps_fire::i, and AlCaHLTBitMon_ParallelJobs::p.

93  {
94  layerClustersPerLayer_.resize(2 * maxlayer_ + 2);
95  // assign all hits in each layer to a cluster core or halo
96  tbb::this_task_arena::isolate([&] {
97  tbb::parallel_for(size_t(0), size_t(2 * maxlayer_ + 2), [&](size_t i) {
98  KDTreeBox bounds(minpos_[i][0], maxpos_[i][0], minpos_[i][1], maxpos_[i][1]);
99  KDTree hit_kdtree;
100  hit_kdtree.build(points_[i], bounds);
101 
102  unsigned int actualLayer =
103  i > maxlayer_ ? (i - (maxlayer_ + 1)) : i; // maps back from index used for KD trees to actual layer
104 
105  double maxdensity = calculateLocalDensity(points_[i], hit_kdtree, actualLayer); // also stores rho (energy
106  // density) for each point (node)
107  // calculate distance to nearest point with higher density storing
108  // distance (delta) and point's index
110  findAndAssignClusters(points_[i], hit_kdtree, maxdensity, bounds, actualLayer, layerClustersPerLayer_[i]);
111  });
112  });
113  //Now that we have the density per point we can store it
114  for (auto const &p : points_) {
115  setDensity(p);
116  }
117 }
KDTreeLinkerAlgo< Hexel > KDTree
std::vector< std::array< float, 2 > > minpos_
std::vector< std::array< float, 2 > > maxpos_
std::vector< std::vector< KDNode > > points_
double calculateLocalDensity(std::vector< KDNode > &, KDTree &, const unsigned int) const
void setDensity(const std::vector< KDNode > &nd)
int findAndAssignClusters(std::vector< KDNode > &, KDTree &, double, KDTreeBox< 2 > &, const unsigned int, std::vector< std::vector< KDNode >> &) const
std::vector< std::vector< std::vector< KDNode > > > layerClustersPerLayer_
double calculateDistanceToHigher(std::vector< KDNode > &) const

◆ populate()

void HGCalImagingAlgo::populate ( const HGCRecHitCollection hits)
overridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 26 of file HGCalImagingAlgo.cc.

References CaloRecHit::detid(), CaloRecHit::energy(), hfClusterShapes_cfi::hits, mps_fire::i, createfilelist::int, SiStripPI::max, SiStripPI::min, position, and Calorimetry_cff::thickness.

26  {
27  // loop over all hits and create the Hexel structure, skip energies below ecut
28 
29  if (dependSensor_) {
30  // for each layer and wafer calculate the thresholds (sigmaNoise and energy)
31  // once
33  }
34 
35  std::vector<bool> firstHit(2 * (maxlayer_ + 1), true);
36  for (unsigned int i = 0; i < hits.size(); ++i) {
37  const HGCRecHit &hgrh = hits[i];
38  DetId detid = hgrh.detid();
39  unsigned int layer = rhtools_.getLayerWithOffset(detid);
40  float thickness = 0.f;
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_) {
46  int thickness_index = rhtools_.getSiThickIndex(detid);
47  if (thickness_index == -1)
48  thickness_index = 3;
49  double storedThreshold = thresholds_[layer - 1][thickness_index];
50  sigmaNoise = sigmaNoise_[layer - 1][thickness_index];
51 
52  if (hgrh.energy() < storedThreshold)
53  continue; // this sets the ZS threshold at ecut times the sigma noise
54  // for the sensor
55  }
56  if (!dependSensor_ && hgrh.energy() < ecut_)
57  continue;
58 
59  // map layers from positive endcap (z) to layer + maxlayer+1 to prevent
60  // mixing up hits from different sides
61  layer += int(rhtools_.zside(detid) > 0) * (maxlayer_ + 1);
62 
63  // determine whether this is a half-hexagon
64  bool isHalf = rhtools_.isHalfCell(detid);
66 
67  // here's were the KDNode is passed its dims arguments - note that these are
68  // *copied* from the Hexel
69  points_[layer].emplace_back(
70  Hexel(hgrh, detid, isHalf, sigmaNoise, thickness, &rhtools_), position.x(), position.y());
71 
72  // for each layer, store the minimum and maximum x and y coordinates for the
73  // KDTreeBox boundaries
74  if (firstHit[layer]) {
75  minpos_[layer][0] = position.x();
76  minpos_[layer][1] = position.y();
77  maxpos_[layer][0] = position.x();
78  maxpos_[layer][1] = position.y();
79  firstHit[layer] = false;
80  } else {
81  minpos_[layer][0] = std::min((float)position.x(), minpos_[layer][0]);
82  minpos_[layer][1] = std::min((float)position.y(), minpos_[layer][1]);
83  maxpos_[layer][0] = std::max((float)position.x(), maxpos_[layer][0]);
84  maxpos_[layer][1] = std::max((float)position.y(), maxpos_[layer][1]);
85  }
86 
87  } // end loop hits
88 }
bool isHalfCell(const DetId &) const
Definition: RecHitTools.cc:415
constexpr const DetId & detid() const
Definition: CaloRecHit.h:33
std::vector< std::vector< double > > thresholds_
std::vector< std::array< float, 2 > > minpos_
int zside(const DetId &id) const
Definition: RecHitTools.cc:163
constexpr float energy() const
Definition: CaloRecHit.h:29
constexpr std::array< uint8_t, layerIndexSize > layer
std::vector< std::array< float, 2 > > maxpos_
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
std::vector< std::vector< KDNode > > points_
std::vector< std::vector< double > > sigmaNoise_
Definition: DetId.h:17
static int position[264][3]
Definition: ReadPGInfo.cc:289
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:205
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:179
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:365

◆ reset()

void HGCalImagingAlgo::reset ( void  )
inlineoverridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 67 of file HGCalImagingAlgo.h.

References clusters_v_, mps_fire::i, layerClustersPerLayer_, maxpos_, minpos_, points_, and edm::swap().

67  {
68  clusters_v_.clear();
69  clusters_v_.shrink_to_fit();
70  layerClustersPerLayer_.clear();
71  layerClustersPerLayer_.shrink_to_fit();
72  for (auto &it : points_) {
73  it.clear();
74  it.shrink_to_fit();
75  std::vector<KDNode>().swap(it);
76  }
77  for (unsigned int i = 0; i < minpos_.size(); i++) {
78  minpos_[i][0] = 0.;
79  minpos_[i][1] = 0.;
80  maxpos_[i][0] = 0.;
81  maxpos_[i][1] = 0.;
82  }
83  }
std::vector< std::array< float, 2 > > minpos_
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
std::vector< std::array< float, 2 > > maxpos_
std::vector< reco::BasicCluster > clusters_v_
std::vector< std::vector< KDNode > > points_
std::vector< std::vector< std::vector< KDNode > > > layerClustersPerLayer_

◆ setDensity()

void HGCalImagingAlgo::setDensity ( const std::vector< KDNode > &  nd)
private

Definition at line 660 of file HGCalImagingAlgo.cc.

References mps_fire::i.

660  {
661  // for each node calculate local density rho and store it
662  for (auto &i : nd) {
663  density_[i.data.detid] = i.data.rho;
664  } // end loop nodes
665 }

◆ shareEnergy()

void HGCalImagingAlgo::shareEnergy ( const std::vector< KDNode > &  incluster,
const std::vector< unsigned > &  seeds,
std::vector< std::vector< double >> &  outclusters 
)
private

Definition at line 531 of file HGCalImagingAlgo.cc.

References data, change_name::diff, MillePedeFileConverter_cfg::e, JetChargeProducer_cfi::exp, DivergingColor::frac, HLT_2022v15_cff::fraction, mps_fire::i, dqmiolumiharvest::j, SiStripPI::max, HLT_2022v15_cff::minFracTot, eostools::move(), perp2(), funct::pow(), DetachedQuadStep_cff::seeds, mathSSE::sqrt(), HLT_2022v15_cff::stoppingTolerance, HGCalImagingAlgo::Hexel::x, HGCalImagingAlgo::Hexel::y, and HGCalImagingAlgo::Hexel::z.

533  {
534  std::vector<bool> isaseed(incluster.size(), false);
535  outclusters.clear();
536  outclusters.resize(seeds.size());
537  std::vector<Point> centroids(seeds.size());
538  std::vector<double> energies(seeds.size());
539 
540  if (seeds.size() == 1) { // short circuit the case of a lone cluster
541  outclusters[0].clear();
542  outclusters[0].resize(incluster.size(), 1.0);
543  return;
544  }
545 
546  // saving seeds
547 
548  // create quick seed lookup
549  for (unsigned i = 0; i < seeds.size(); ++i) {
550  isaseed[seeds[i]] = true;
551  }
552 
553  // initialize clusters to be shared
554  // centroids start off at seed positions
555  // seeds always have fraction 1.0, to stabilize fit
556  // initializing fit
557  for (unsigned i = 0; i < seeds.size(); ++i) {
558  outclusters[i].resize(incluster.size(), 0.0);
559  for (unsigned j = 0; j < incluster.size(); ++j) {
560  if (j == seeds[i]) {
561  outclusters[i][j] = 1.0;
562  centroids[i] = math::XYZPoint(incluster[j].data.x, incluster[j].data.y, incluster[j].data.z);
563  energies[i] = incluster[j].data.weight;
564  }
565  }
566  }
567 
568  // run the fit while we are less than max iterations, and clusters are still
569  // moving
570  const double minFracTot = 1e-20;
571  unsigned iter = 0;
572  const unsigned iterMax = 50;
574  const double stoppingTolerance = 1e-8;
575  const auto numberOfSeeds = seeds.size();
576  auto toleranceScaling = numberOfSeeds > 2 ? (numberOfSeeds - 1) * (numberOfSeeds - 1) : 1;
577  std::vector<Point> prevCentroids;
578  std::vector<double> frac(numberOfSeeds), dist2(numberOfSeeds);
579  while (iter++ < iterMax && diff > stoppingTolerance * toleranceScaling) {
580  for (unsigned i = 0; i < incluster.size(); ++i) {
581  const Hexel &ihit = incluster[i].data;
582  double fracTot(0.0);
583  for (unsigned j = 0; j < numberOfSeeds; ++j) {
584  double fraction = 0.0;
585  double d2 = (std::pow(ihit.x - centroids[j].x(), 2.0) + std::pow(ihit.y - centroids[j].y(), 2.0) +
586  std::pow(ihit.z - centroids[j].z(), 2.0)) /
587  sigma2_;
588  dist2[j] = d2;
589  // now we set the fractions up based on hit type
590  if (i == seeds[j]) { // this cluster's seed
591  fraction = 1.0;
592  } else if (isaseed[i]) {
593  fraction = 0.0;
594  } else {
595  fraction = energies[j] * std::exp(-0.5 * d2);
596  }
597  fracTot += fraction;
598  frac[j] = fraction;
599  }
600  // now that we have calculated all fractions for all hits
601  // assign the new fractions
602  for (unsigned j = 0; j < numberOfSeeds; ++j) {
603  if (fracTot > minFracTot || (i == seeds[j] && fracTot > 0.0)) {
604  outclusters[j][i] = frac[j] / fracTot;
605  } else {
606  outclusters[j][i] = 0.0;
607  }
608  }
609  }
610 
611  // save previous centroids
612  prevCentroids = std::move(centroids);
613  // finally update the position of the centroids from the last iteration
614  centroids.resize(numberOfSeeds);
615  double diff2 = 0.0;
616  for (unsigned i = 0; i < numberOfSeeds; ++i) {
617  centroids[i] = calculatePositionWithFraction(incluster, outclusters[i]);
618  energies[i] = calculateEnergyWithFraction(incluster, outclusters[i]);
619  // calculate convergence parameters
620  const double delta2 = (prevCentroids[i] - centroids[i]).perp2();
621  diff2 = std::max(delta2, diff2);
622  }
623  // update convergance parameter outside loop
624  diff = std::sqrt(diff2);
625  }
626 }
double calculateEnergyWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
math::XYZPoint calculatePositionWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
T perp2() const
Squared magnitude of transverse component.
T sqrt(T t)
Definition: SSEVec.h:19
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511

◆ sort_by_delta()

std::vector<size_t> HGCalImagingAlgo::sort_by_delta ( const std::vector< KDNode > &  v) const
inlineprivate

Definition at line 225 of file HGCalImagingAlgo.h.

References testProducerWithPsetDescEmpty_cfi::i1, testProducerWithPsetDescEmpty_cfi::i2, heavyIonCSV_trainingSettings::idx, jetsAK4_CHS_cff::sort, and findQualityFiles::v.

225  {
226  std::vector<size_t> idx(v.size());
227  std::iota(std::begin(idx), std::end(idx), 0);
228  sort(idx.begin(), idx.end(), [&v](size_t i1, size_t i2) { return v[i1].data.delta > v[i2].data.delta; });
229  return idx;
230  }

Member Data Documentation

◆ clusters_v_

std::vector<reco::BasicCluster> HGCalImagingAlgo::clusters_v_
private

Definition at line 138 of file HGCalImagingAlgo.h.

Referenced by reset().

◆ dEdXweights_

std::vector<double> HGCalImagingAlgo::dEdXweights_
private

Definition at line 145 of file HGCalImagingAlgo.h.

◆ density_

Density HGCalImagingAlgo::density_
private

Definition at line 141 of file HGCalImagingAlgo.h.

◆ dependSensor_

bool HGCalImagingAlgo::dependSensor_
private

Definition at line 144 of file HGCalImagingAlgo.h.

◆ ecut_

double HGCalImagingAlgo::ecut_
private

Definition at line 132 of file HGCalImagingAlgo.h.

◆ fcPerEle_

double HGCalImagingAlgo::fcPerEle_
private

Definition at line 148 of file HGCalImagingAlgo.h.

◆ fcPerMip_

std::vector<double> HGCalImagingAlgo::fcPerMip_
private

Definition at line 147 of file HGCalImagingAlgo.h.

◆ initialized_

bool HGCalImagingAlgo::initialized_
private

Definition at line 155 of file HGCalImagingAlgo.h.

◆ kappa_

double HGCalImagingAlgo::kappa_
private

Definition at line 129 of file HGCalImagingAlgo.h.

◆ layerClustersPerLayer_

std::vector<std::vector<std::vector<KDNode> > > HGCalImagingAlgo::layerClustersPerLayer_
private

Definition at line 223 of file HGCalImagingAlgo.h.

Referenced by reset().

◆ maxpos_

std::vector<std::array<float, 2> > HGCalImagingAlgo::maxpos_
private

Definition at line 236 of file HGCalImagingAlgo.h.

Referenced by reset().

◆ minpos_

std::vector<std::array<float, 2> > HGCalImagingAlgo::minpos_
private

Definition at line 235 of file HGCalImagingAlgo.h.

Referenced by reset().

◆ noiseMip_

double HGCalImagingAlgo::noiseMip_
private

Definition at line 150 of file HGCalImagingAlgo.h.

◆ nonAgedNoises_

std::vector<double> HGCalImagingAlgo::nonAgedNoises_
private

Definition at line 149 of file HGCalImagingAlgo.h.

◆ points_

std::vector<std::vector<KDNode> > HGCalImagingAlgo::points_
private

Definition at line 232 of file HGCalImagingAlgo.h.

Referenced by reset().

◆ positionDeltaRho_c_

std::vector<double> HGCalImagingAlgo::positionDeltaRho_c_
private

Definition at line 125 of file HGCalImagingAlgo.h.

◆ sigma2_

double HGCalImagingAlgo::sigma2_
private

Definition at line 135 of file HGCalImagingAlgo.h.

◆ sigmaNoise_

std::vector<std::vector<double> > HGCalImagingAlgo::sigmaNoise_
private

Definition at line 152 of file HGCalImagingAlgo.h.

◆ thicknessCorrection_

std::vector<double> HGCalImagingAlgo::thicknessCorrection_
private

Definition at line 146 of file HGCalImagingAlgo.h.

◆ thresholds_

std::vector<std::vector<double> > HGCalImagingAlgo::thresholds_
private

Definition at line 151 of file HGCalImagingAlgo.h.

◆ thresholdW0_

std::vector<double> HGCalImagingAlgo::thresholdW0_
private

Definition at line 124 of file HGCalImagingAlgo.h.

◆ vecDeltas_

std::vector<double> HGCalImagingAlgo::vecDeltas_
private

Definition at line 128 of file HGCalImagingAlgo.h.