CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
ticl::PatternRecognitionbyCLUE3D< TILES > Class Template Referencefinal

#include <PatternRecognitionbyCLUE3D.h>

Inheritance diagram for ticl::PatternRecognitionbyCLUE3D< TILES >:
ticl::PatternRecognitionAlgoBaseT< TILES >

Classes

struct  ClustersOnLayer
 

Public Member Functions

void energyRegressionAndID (const std::vector< reco::CaloCluster > &layerClusters, const tensorflow::Session *, std::vector< Trackster > &result)
 
void filter (std::vector< Trackster > &output, const std::vector< Trackster > &inTracksters, const typename PatternRecognitionAlgoBaseT< TILES >::Inputs &input, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation) override
 
void makeTracksters (const typename PatternRecognitionAlgoBaseT< TILES >::Inputs &input, std::vector< Trackster > &result, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation) override
 
 PatternRecognitionbyCLUE3D (const edm::ParameterSet &conf, edm::ConsumesCollector)
 
 ~PatternRecognitionbyCLUE3D () override=default
 
- Public Member Functions inherited from ticl::PatternRecognitionAlgoBaseT< TILES >
virtual void filter (std::vector< Trackster > &output, const std::vector< Trackster > &inTracksters, const Inputs &input, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation)=0
 
virtual void makeTracksters (const Inputs &input, std::vector< Trackster > &result, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation)=0
 
 PatternRecognitionAlgoBaseT (const edm::ParameterSet &conf, edm::ConsumesCollector)
 
virtual ~PatternRecognitionAlgoBaseT ()
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &iDesc)
 

Private Member Functions

void calculateDistanceToHigher (const TILES &, const int layerId, const std::vector< std::pair< int, int >> &)
 
void calculateLocalDensity (const TILES &, const int layerId, const std::vector< std::pair< int, int >> &)
 
void dumpClusters (const TILES &tiles, const std::vector< std::pair< int, int >> &layerIdx2layerandSoa, const int) const
 
void dumpTiles (const TILES &) const
 
void dumpTracksters (const std::vector< std::pair< int, int >> &layerIdx2layerandSoa, const int, const std::vector< Trackster > &) const
 
int findAndAssignTracksters (const TILES &, const std::vector< std::pair< int, int >> &)
 
void reset ()
 

Private Attributes

edm::ESGetToken< CaloGeometry, CaloGeometryRecordcaloGeomToken_
 
std::vector< ClustersOnLayerclusters_
 
const bool computeLocalTime_
 
const std::vector< double > criticalDensity_
 
const std::vector< double > criticalEtaPhiDistance_
 
const std::vector< double > criticalSelfDensity_
 
const std::vector< double > criticalXYDistance_
 
const std::vector< int > criticalZDistanceLyr_
 
const float cutHadProb_
 
const std::vector< double > densityEtaPhiDistanceSqr_
 
const bool densityOnSameLayer_
 
const std::vector< int > densitySiblingLayers_
 
const std::vector< double > densityXYDistanceSqr_
 
const bool doPidCut_
 
const std::vector< int > filter_on_categories_
 
const std::vector< double > kernelDensityFactor_
 
std::vector< float > layersPosZ_
 
const std::vector< int > minNumLayerCluster_
 
const bool nearestHigherOnSameLayer_
 
const std::vector< double > outlierMultiplier_
 
const bool rescaleDensityByZ_
 
hgcal::RecHitTools rhtools_
 
std::vector< int > tracksterSeedAlgoId_
 
const bool useAbsoluteProjectiveScale_
 
const bool useClusterDimensionXY_
 
const bool usePCACleaning_
 

Additional Inherited Members

- Protected Attributes inherited from ticl::PatternRecognitionAlgoBaseT< TILES >
int algo_verbosity_
 

Detailed Description

template<typename TILES>
class ticl::PatternRecognitionbyCLUE3D< TILES >

Definition at line 12 of file PatternRecognitionbyCLUE3D.h.

Constructor & Destructor Documentation

◆ PatternRecognitionbyCLUE3D()

template<typename TILES >
PatternRecognitionbyCLUE3D::PatternRecognitionbyCLUE3D ( const edm::ParameterSet conf,
edm::ConsumesCollector  iC 
)

Definition at line 21 of file PatternRecognitionbyCLUE3D.cc.

24  criticalDensity_(conf.getParameter<std::vector<double>>("criticalDensity")),
25  criticalSelfDensity_(conf.getParameter<std::vector<double>>("criticalSelfDensity")),
26  densitySiblingLayers_(conf.getParameter<std::vector<int>>("densitySiblingLayers")),
27  densityEtaPhiDistanceSqr_(conf.getParameter<std::vector<double>>("densityEtaPhiDistanceSqr")),
28  densityXYDistanceSqr_(conf.getParameter<std::vector<double>>("densityXYDistanceSqr")),
29  kernelDensityFactor_(conf.getParameter<std::vector<double>>("kernelDensityFactor")),
30  densityOnSameLayer_(conf.getParameter<bool>("densityOnSameLayer")),
31  nearestHigherOnSameLayer_(conf.getParameter<bool>("nearestHigherOnSameLayer")),
32  useAbsoluteProjectiveScale_(conf.getParameter<bool>("useAbsoluteProjectiveScale")),
33  useClusterDimensionXY_(conf.getParameter<bool>("useClusterDimensionXY")),
34  rescaleDensityByZ_(conf.getParameter<bool>("rescaleDensityByZ")),
35  criticalEtaPhiDistance_(conf.getParameter<std::vector<double>>("criticalEtaPhiDistance")),
36  criticalXYDistance_(conf.getParameter<std::vector<double>>("criticalXYDistance")),
37  criticalZDistanceLyr_(conf.getParameter<std::vector<int>>("criticalZDistanceLyr")),
38  outlierMultiplier_(conf.getParameter<std::vector<double>>("outlierMultiplier")),
39  minNumLayerCluster_(conf.getParameter<std::vector<int>>("minNumLayerCluster")),
40  doPidCut_(conf.getParameter<bool>("doPidCut")),
41  cutHadProb_(conf.getParameter<double>("cutHadProb")),
42  computeLocalTime_(conf.getParameter<bool>("computeLocalTime")),
43  usePCACleaning_(conf.getParameter<bool>("usePCACleaning")){};
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
const std::vector< double > densityEtaPhiDistanceSqr_
const std::vector< double > densityXYDistanceSqr_
const std::vector< double > outlierMultiplier_
const std::vector< int > densitySiblingLayers_
const std::vector< double > kernelDensityFactor_
const std::vector< double > criticalDensity_
const std::vector< double > criticalXYDistance_
const std::vector< double > criticalEtaPhiDistance_
const std::vector< int > criticalZDistanceLyr_
const std::vector< double > criticalSelfDensity_

◆ ~PatternRecognitionbyCLUE3D()

template<typename TILES >
ticl::PatternRecognitionbyCLUE3D< TILES >::~PatternRecognitionbyCLUE3D ( )
overridedefault

Member Function Documentation

◆ calculateDistanceToHigher()

template<typename TILES >
void PatternRecognitionbyCLUE3D::calculateDistanceToHigher ( const TILES &  tiles,
const int  layerId,
const std::vector< std::pair< int, int >> &  layerIdx2layerandSoa 
)
private

Definition at line 547 of file PatternRecognitionbyCLUE3D.cc.

References funct::abs(), ticl::Advanced, hltPFPuppi_cfi::algoId, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), HLT_2024v14_cff::delta_phi, reco::deltaPhi(), reco::deltaR2(), mps_fire::i, hcalRecHitTable_cff::ieta, hcalRecHitTable_cff::iphi, WZElectronSkims53X_cff::max, allConversions_cfi::maxDelta, SiStripPI::min, HLT_2024v14_cff::minLayer, L1TMuonDQMOffline_cfi::nEtaBins, ecaldqm::binning::nPhiBins, electrons_cff::numberOfClusters, HLT_IsoTrack_cff::offset, trackingPOGFilters_cfi::phiWindow, and mathSSE::sqrt().

548  {
551  auto &clustersOnLayer = clusters_[layerId];
552  unsigned int numberOfClusters = clustersOnLayer.x.size();
553 
554  auto distanceSqr = [](float r0, float r1, float phi0, float phi1) -> float {
555  auto delta_phi = reco::deltaPhi(phi0, phi1);
556  return (r0 - r1) * (r0 - r1) + r1 * r1 * delta_phi * delta_phi;
557  };
558 
559  for (unsigned int i = 0; i < numberOfClusters; i++) {
561  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
562  << "Starting searching nearestHigher on " << layerId << " with rho: " << clustersOnLayer.rho[i]
563  << " at eta, phi: " << tiles[layerId].etaBin(clustersOnLayer.eta[i]) << ", "
564  << tiles[layerId].phiBin(clustersOnLayer.phi[i]);
565  }
566  // We need to partition the two sides of the HGCAL detector
567  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
568  int minLayer = 0;
569  int maxLayer = 2 * lastLayerPerSide - 1;
570  auto algoId = clustersOnLayer.algoId[i];
571  if (layerId < lastLayerPerSide) {
573  maxLayer = std::min(layerId + densitySiblingLayers_[algoId], lastLayerPerSide - 1);
574  } else {
575  minLayer = std::max(layerId - densitySiblingLayers_[algoId], lastLayerPerSide + 1);
576  maxLayer = std::min(layerId + densitySiblingLayers_[algoId], maxLayer);
577  }
579  float i_delta = maxDelta;
580  std::pair<int, int> i_nearestHigher(-1, -1);
581  std::pair<float, int> nearest_distances(maxDelta, std::numeric_limits<int>::max());
582  for (int currentLayer = minLayer; currentLayer <= maxLayer; currentLayer++) {
583  if (!nearestHigherOnSameLayer_ && (layerId == currentLayer))
584  continue;
585  const auto &tileOnLayer = tiles[currentLayer];
586  int etaWindow = 1;
587  int phiWindow = 1;
588  int etaBinMin = std::max(tileOnLayer.etaBin(clustersOnLayer.eta[i]) - etaWindow, 0);
589  int etaBinMax = std::min(tileOnLayer.etaBin(clustersOnLayer.eta[i]) + etaWindow, nEtaBin - 1);
590  int phiBinMin = tileOnLayer.phiBin(clustersOnLayer.phi[i]) - phiWindow;
591  int phiBinMax = tileOnLayer.phiBin(clustersOnLayer.phi[i]) + phiWindow;
592  for (int ieta = etaBinMin; ieta <= etaBinMax; ++ieta) {
593  auto offset = ieta * nPhiBin;
594  for (int iphi_it = phiBinMin; iphi_it <= phiBinMax; ++iphi_it) {
595  int iphi = ((iphi_it % nPhiBin + nPhiBin) % nPhiBin);
597  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
598  << "Searching nearestHigher on " << currentLayer << " eta, phi: " << ieta << ", " << iphi_it << " "
599  << iphi << " " << offset << " " << (offset + iphi);
600  }
601  for (auto otherClusterIdx : tileOnLayer[offset + iphi]) {
602  auto const &layerandSoa = layerIdx2layerandSoa[otherClusterIdx];
603  // Skip masked layer clusters
604  if ((layerandSoa.first == -1) && (layerandSoa.second == -1))
605  continue;
606  auto const &clustersOnOtherLayer = clusters_[layerandSoa.first];
607  auto dist = maxDelta;
608  auto dist_transverse = maxDelta;
609  int dist_layers = std::abs(layerandSoa.first - layerId);
611  dist_transverse = distanceSqr(clustersOnLayer.r_over_absz[i] * clustersOnLayer.z[i],
612  clustersOnOtherLayer.r_over_absz[layerandSoa.second] * clustersOnLayer.z[i],
613  clustersOnLayer.phi[i],
614  clustersOnOtherLayer.phi[layerandSoa.second]);
615  // Add Z-scale to the final distance
616  dist = dist_transverse;
617  } else {
618  dist = reco::deltaR2(clustersOnLayer.eta[i],
619  clustersOnLayer.phi[i],
620  clustersOnOtherLayer.eta[layerandSoa.second],
621  clustersOnOtherLayer.phi[layerandSoa.second]);
622  dist_transverse = dist;
623  }
624  bool foundHigher = (clustersOnOtherLayer.rho[layerandSoa.second] > clustersOnLayer.rho[i]) ||
625  (clustersOnOtherLayer.rho[layerandSoa.second] == clustersOnLayer.rho[i] &&
626  clustersOnOtherLayer.layerClusterOriginalIdx[layerandSoa.second] >
627  clustersOnLayer.layerClusterOriginalIdx[i]);
629  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
630  << "Searching nearestHigher on " << currentLayer
631  << " with rho: " << clustersOnOtherLayer.rho[layerandSoa.second]
632  << " on layerIdxInSOA: " << layerandSoa.first << ", " << layerandSoa.second
633  << " with distance: " << sqrt(dist) << " foundHigher: " << foundHigher;
634  }
635  if (foundHigher && dist <= i_delta) {
636  // update i_delta
637  i_delta = dist;
638  nearest_distances = std::make_pair(sqrt(dist_transverse), dist_layers);
639  // update i_nearestHigher
640  i_nearestHigher = layerandSoa;
641  }
642  } // End of loop on clusters
643  } // End of loop on phi bins
644  } // End of loop on eta bins
645  } // End of loop on layers
646 
647  clustersOnLayer.delta[i] = nearest_distances;
648  clustersOnLayer.nearestHigher[i] = i_nearestHigher;
649  } // End of loop on clusters
650 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_
const std::vector< int > densitySiblingLayers_
T sqrt(T t)
Definition: SSEVec.h:23
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static constexpr int nPhiBins
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:80

◆ calculateLocalDensity()

template<typename TILES >
void PatternRecognitionbyCLUE3D::calculateLocalDensity ( const TILES &  tiles,
const int  layerId,
const std::vector< std::pair< int, int >> &  layerIdx2layerandSoa 
)
private

Definition at line 381 of file PatternRecognitionbyCLUE3D.cc.

References funct::abs(), ticl::Advanced, hltPFPuppi_cfi::algoId, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), HLT_2024v14_cff::delta_phi, reco::deltaPhi(), reco::deltaR2(), f, mps_fire::i, hcalRecHitTable_cff::ieta, ALPAKA_ACCELERATOR_NAMESPACE::caPixelDoublets::if(), hcalRecHitTable_cff::iphi, WZElectronSkims53X_cff::max, SiStripPI::min, HLT_2024v14_cff::minLayer, L1TMuonDQMOffline_cfi::nEtaBins, ecaldqm::binning::nPhiBins, electrons_cff::numberOfClusters, HLT_IsoTrack_cff::offset, trackingPOGFilters_cfi::phiWindow, mathSSE::sqrt(), testProducerWithPsetDescEmpty_cfi::y1, and testProducerWithPsetDescEmpty_cfi::y2.

382  {
385  auto &clustersOnLayer = clusters_[layerId];
386  unsigned int numberOfClusters = clustersOnLayer.x.size();
387 
388  auto isReachable = [](float r0, float r1, float phi0, float phi1, float delta_sqr) -> bool {
389  auto delta_phi = reco::deltaPhi(phi0, phi1);
390  return (r0 - r1) * (r0 - r1) + r1 * r1 * delta_phi * delta_phi < delta_sqr;
391  };
392  auto distance_debug = [&](float x1, float x2, float y1, float y2) -> float {
393  return sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2));
394  };
395 
396  for (unsigned int i = 0; i < numberOfClusters; i++) {
397  auto algoId = clustersOnLayer.algoId[i];
398  // We need to partition the two sides of the HGCAL detector
399  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
400  int minLayer = 0;
401  int maxLayer = 2 * lastLayerPerSide - 1;
402  if (layerId < lastLayerPerSide) {
404  maxLayer = std::min(layerId + densitySiblingLayers_[algoId], lastLayerPerSide - 1);
405  } else {
406  minLayer = std::max(layerId - densitySiblingLayers_[algoId], lastLayerPerSide);
407  maxLayer = std::min(layerId + densitySiblingLayers_[algoId], maxLayer);
408  }
409  float deltaLayersZ = std::abs(layersPosZ_[maxLayer % lastLayerPerSide] - layersPosZ_[minLayer % lastLayerPerSide]);
410 
411  for (int currentLayer = minLayer; currentLayer <= maxLayer; currentLayer++) {
413  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "RefLayer: " << layerId << " SoaIDX: " << i;
414  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "NextLayer: " << currentLayer;
415  }
416  const auto &tileOnLayer = tiles[currentLayer];
417  bool onSameLayer = (currentLayer == layerId);
419  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "onSameLayer: " << onSameLayer;
420  }
421  const int etaWindow = 2;
422  const int phiWindow = 2;
423  int etaBinMin = std::max(tileOnLayer.etaBin(clustersOnLayer.eta[i]) - etaWindow, 0);
424  int etaBinMax = std::min(tileOnLayer.etaBin(clustersOnLayer.eta[i]) + etaWindow, nEtaBin - 1);
425  int phiBinMin = tileOnLayer.phiBin(clustersOnLayer.phi[i]) - phiWindow;
426  int phiBinMax = tileOnLayer.phiBin(clustersOnLayer.phi[i]) + phiWindow;
428  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "eta: " << clustersOnLayer.eta[i];
429  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "phi: " << clustersOnLayer.phi[i];
430  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "etaBinMin: " << etaBinMin << ", etaBinMax: " << etaBinMax;
431  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "phiBinMin: " << phiBinMin << ", phiBinMax: " << phiBinMax;
432  }
433  for (int ieta = etaBinMin; ieta <= etaBinMax; ++ieta) {
434  auto offset = ieta * nPhiBin;
436  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "offset: " << offset;
437  }
438  for (int iphi_it = phiBinMin; iphi_it <= phiBinMax; ++iphi_it) {
439  int iphi = ((iphi_it % nPhiBin + nPhiBin) % nPhiBin);
441  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "iphi: " << iphi;
442  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
443  << "Entries in tileBin: " << tileOnLayer[offset + iphi].size();
444  }
445  for (auto otherClusterIdx : tileOnLayer[offset + iphi]) {
446  auto const &layerandSoa = layerIdx2layerandSoa[otherClusterIdx];
447  // Skip masked layer clusters
448  if ((layerandSoa.first == -1) && (layerandSoa.second == -1)) {
450  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Skipping masked layerIdx " << otherClusterIdx;
451  }
452  continue;
453  }
454  auto const &clustersLayer = clusters_[layerandSoa.first];
456  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
457  << "OtherLayer: " << layerandSoa.first << " SoaIDX: " << layerandSoa.second;
458  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "OtherEta: " << clustersLayer.eta[layerandSoa.second];
459  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "OtherPhi: " << clustersLayer.phi[layerandSoa.second];
460  }
461 
462  bool onSameCluster = clustersOnLayer.layerClusterOriginalIdx[i] == otherClusterIdx;
463  if (onSameLayer && !densityOnSameLayer_ && !onSameCluster) {
465  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
466  << "Skipping different cluster " << otherClusterIdx << "in the same layer " << currentLayer;
467  }
468  continue;
469  }
470 
471  bool reachable = false;
474  reachable = isReachable(clustersOnLayer.r_over_absz[i] * clustersOnLayer.z[i],
475  clustersLayer.r_over_absz[layerandSoa.second] * clustersOnLayer.z[i],
476  clustersOnLayer.phi[i],
477  clustersLayer.phi[layerandSoa.second],
478  clustersOnLayer.radius[i] * clustersOnLayer.radius[i]);
479  } else {
480  // Still differentiate between silicon and Scintillator.
481  // Scintillator has yet to be studied further.
482  if (clustersOnLayer.isSilicon[i]) {
483  reachable = isReachable(clustersOnLayer.r_over_absz[i] * clustersOnLayer.z[i],
484  clustersLayer.r_over_absz[layerandSoa.second] * clustersOnLayer.z[i],
485  clustersOnLayer.phi[i],
486  clustersLayer.phi[layerandSoa.second],
488  } else {
489  reachable = isReachable(clustersOnLayer.r_over_absz[i] * clustersOnLayer.z[i],
490  clustersLayer.r_over_absz[layerandSoa.second] * clustersOnLayer.z[i],
491  clustersOnLayer.phi[i],
492  clustersLayer.phi[layerandSoa.second],
493  clustersOnLayer.radius[i] * clustersOnLayer.radius[i]);
494  }
495  }
496  } else {
497  reachable = (reco::deltaR2(clustersOnLayer.eta[i],
498  clustersOnLayer.phi[i],
499  clustersLayer.eta[layerandSoa.second],
500  clustersLayer.phi[layerandSoa.second]) < densityEtaPhiDistanceSqr_[algoId]);
501  }
503  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Distance[eta,phi]: "
504  << reco::deltaR2(clustersOnLayer.eta[i],
505  clustersOnLayer.phi[i],
506  clustersLayer.eta[layerandSoa.second],
507  clustersLayer.phi[layerandSoa.second]);
508  auto dist = distance_debug(
509  clustersOnLayer.r_over_absz[i],
510  clustersLayer.r_over_absz[layerandSoa.second],
511  clustersOnLayer.r_over_absz[i] * std::abs(clustersOnLayer.phi[i]),
512  clustersLayer.r_over_absz[layerandSoa.second] * std::abs(clustersLayer.phi[layerandSoa.second]));
513  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Distance[cm]: " << (dist * clustersOnLayer.z[i]);
514  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
515  << "Energy Other: " << clustersLayer.energy[layerandSoa.second];
516  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Cluster radius: " << clustersOnLayer.radius[i];
517  }
518  if (reachable) {
519  float factor_same_layer_different_cluster = (onSameLayer && !densityOnSameLayer_) ? 0.f : 1.f;
520  auto energyToAdd = (clustersOnLayer.layerClusterOriginalIdx[i] == otherClusterIdx
521  ? 1.f
522  : kernelDensityFactor_[algoId] * factor_same_layer_different_cluster) *
523  clustersLayer.energy[layerandSoa.second];
524  clustersOnLayer.rho[i] += energyToAdd;
525  clustersOnLayer.z_extension[i] = deltaLayersZ;
527  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
528  << "Adding " << energyToAdd << " partial " << clustersOnLayer.rho[i];
529  }
530  }
531  } // end of loop on possible compatible clusters
532  } // end of loop over phi-bin region
533  } // end of loop over eta-bin region
534  } // end of loop on the sibling layers
535  if (rescaleDensityByZ_) {
537  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
538  << "Rescaling original density: " << clustersOnLayer.rho[i] << " by Z: " << deltaLayersZ
539  << " to final density/cm: " << clustersOnLayer.rho[i] / deltaLayersZ;
540  }
541  clustersOnLayer.rho[i] /= deltaLayersZ;
542  }
543  } // end of loop over clusters on this layer
544 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
Log< level::Info, true > LogVerbatim
const std::vector< double > densityEtaPhiDistanceSqr_
const std::vector< double > densityXYDistanceSqr_
std::vector< ClustersOnLayer > clusters_
const std::vector< int > densitySiblingLayers_
T sqrt(T t)
Definition: SSEVec.h:23
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const std::vector< double > kernelDensityFactor_
double f[11][100]
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static constexpr int nPhiBins
if(threadIdxLocalY==0 &&threadIdxLocalX==0)
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:80

◆ dumpClusters()

template<typename TILES >
void PatternRecognitionbyCLUE3D::dumpClusters ( const TILES &  tiles,
const std::vector< std::pair< int, int >> &  layerIdx2layerandSoa,
const int  eventNumber 
) const
private

Definition at line 102 of file PatternRecognitionbyCLUE3D.cc.

References ticl::Advanced, alignBH_cfg::fixed, EgammaValidation_cff::num, and findQualityFiles::v.

104  {
106  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "[evt, lyr, Seed, x, y, z, r/|z|, eta, phi, "
107  "etab, phib, cells, enrgy, e/rho, rho, z_ext, "
108  " dlt_tr, dlt_lyr, "
109  " nestHL, nestHSoaIdx, radius, clIdx, lClOrigIdx, SOAidx";
110  }
111 
112  for (unsigned int layer = 0; layer < clusters_.size(); layer++) {
113  auto const &thisLayer = clusters_[layer];
114  int num = 0;
115  for (auto v : thisLayer.x) {
117  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
118  << std::setw(4) << eventNumber << ", " << std::setw(3) << layer << ", " << std::setw(4)
119  << thisLayer.isSeed[num] << ", " << std::setprecision(3) << std::fixed << v << ", " << thisLayer.y[num]
120  << ", " << thisLayer.z[num] << ", " << thisLayer.r_over_absz[num] << ", " << thisLayer.eta[num] << ", "
121  << thisLayer.phi[num] << ", " << std::setw(5) << tiles[layer].etaBin(thisLayer.eta[num]) << ", "
122  << std::setw(5) << tiles[layer].phiBin(thisLayer.phi[num]) << ", " << std::setw(4) << thisLayer.cells[num]
123  << ", " << std::setprecision(3) << thisLayer.energy[num] << ", "
124  << (thisLayer.energy[num] / thisLayer.rho[num]) << ", " << thisLayer.rho[num] << ", "
125  << thisLayer.z_extension[num] << ", " << std::scientific << thisLayer.delta[num].first << ", "
126  << std::setw(10) << thisLayer.delta[num].second << ", " << std::setw(5)
127  << thisLayer.nearestHigher[num].first << ", " << std::setw(10) << thisLayer.nearestHigher[num].second
128  << ", " << std::defaultfloat << std::setprecision(3) << thisLayer.radius[num] << ", " << std::setw(5)
129  << thisLayer.clusterIndex[num] << ", " << std::setw(4) << thisLayer.layerClusterOriginalIdx[num] << ", "
130  << std::setw(4) << num << ", ClusterInfo";
131  }
132  ++num;
133  }
134  }
135  for (unsigned int lcIdx = 0; lcIdx < layerIdx2layerandSoa.size(); lcIdx++) {
136  auto const &layerandSoa = layerIdx2layerandSoa[lcIdx];
137  // Skip masked layer clusters
138  if ((layerandSoa.first == -1) && (layerandSoa.second == -1))
139  continue;
141  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
142  << "lcIdx: " << lcIdx << " on Layer: " << layerandSoa.first << " SOA: " << layerandSoa.second;
143  }
144  }
145 }
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_

◆ dumpTiles()

template<typename TILES >
void PatternRecognitionbyCLUE3D::dumpTiles ( const TILES &  tiles) const
private

Definition at line 45 of file PatternRecognitionbyCLUE3D.cc.

References ticl::Advanced, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), relativeConstraints::empty, hcalRecHitTable_cff::ieta, hcalRecHitTable_cff::iphi, L1TMuonDQMOffline_cfi::nEtaBins, ecaldqm::binning::nPhiBins, and HLT_IsoTrack_cff::offset.

45  {
48  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
49  int maxLayer = 2 * lastLayerPerSide - 1;
50  for (int layer = 0; layer <= maxLayer; layer++) {
51  for (int ieta = 0; ieta < nEtaBin; ieta++) {
52  auto offset = ieta * nPhiBin;
53  for (int phi = 0; phi < nPhiBin; phi++) {
54  int iphi = ((phi % nPhiBin + nPhiBin) % nPhiBin);
55  if (!tiles[layer][offset + iphi].empty()) {
57  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Layer: " << layer << " ieta: " << ieta << " phi: " << phi
58  << " " << tiles[layer][offset + iphi].size();
59  }
60  }
61  }
62  }
63  }
64 }
Log< level::Info, true > LogVerbatim
static constexpr int nPhiBins
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:80

◆ dumpTracksters()

template<typename TILES >
void PatternRecognitionbyCLUE3D::dumpTracksters ( const std::vector< std::pair< int, int >> &  layerIdx2layerandSoa,
const int  eventNumber,
const std::vector< Trackster > &  tracksters 
) const
private

Definition at line 67 of file PatternRecognitionbyCLUE3D.cc.

References ticl::Advanced, ticl::Trackster::charged_hadron, ticl::Trackster::electron, ticl::Trackster::neutral_hadron, EgammaValidation_cff::num, ticl::Trackster::photon, AlCaHLTBitMon_QueryRunRegistry::string, submitPVValidationJobs::t, HitToTracksterAssociation_cfi::tracksters, and findQualityFiles::v.

69  {
71  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
72  << "[evt, tracksterId, cells, prob_photon, prob_ele, prob_chad, prob_nhad, layer_i, x_i, y_i, eta_i, phi_i, "
73  "energy_i, radius_i, rho_i, z_extension, delta_tr, delta_lyr, isSeed_i";
74  }
75 
76  int num = 0;
77  const std::string sep(", ");
78  for (auto const &t : tracksters) {
79  for (auto v : t.vertices()) {
80  auto [lyrIdx, soaIdx] = layerIdx2layerandSoa[v];
81  auto const &thisLayer = clusters_[lyrIdx];
83  edm::LogVerbatim("PatternRecognitionbyCLUE3D_NTP")
84  << std::setw(4) << eventNumber << sep << std::setw(4) << num << sep << std::setw(4) << t.vertices().size()
85  << sep << std::setw(8) << t.id_probability(ticl::Trackster::ParticleType::photon) << sep << std::setw(8)
86  << t.id_probability(ticl::Trackster::ParticleType::electron) << sep << std::setw(8)
87  << t.id_probability(ticl::Trackster::ParticleType::charged_hadron) << sep << std::setw(8)
88  << t.id_probability(ticl::Trackster::ParticleType::neutral_hadron) << sep << std::setw(4) << lyrIdx << sep
89  << std::setw(10) << thisLayer.x[soaIdx] << sep << std::setw(10) << thisLayer.y[soaIdx] << sep
90  << std::setw(10) << thisLayer.eta[soaIdx] << sep << std::setw(10) << thisLayer.phi[soaIdx] << sep
91  << std::setw(10) << thisLayer.energy[soaIdx] << sep << std::setw(10) << thisLayer.radius[soaIdx] << sep
92  << std::setw(10) << thisLayer.rho[soaIdx] << sep << std::setw(10) << thisLayer.z_extension[soaIdx] << sep
93  << std::setw(10) << thisLayer.delta[soaIdx].first << sep << std::setw(10) << thisLayer.delta[soaIdx].second
94  << sep << std::setw(4) << thisLayer.isSeed[soaIdx];
95  }
96  }
97  num++;
98  }
99 }
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_

◆ energyRegressionAndID()

template<typename TILES >
void ticl::PatternRecognitionbyCLUE3D< TILES >::energyRegressionAndID ( const std::vector< reco::CaloCluster > &  layerClusters,
const tensorflow::Session *  ,
std::vector< Trackster > &  result 
)

◆ fillPSetDescription()

template<typename TILES >
void PatternRecognitionbyCLUE3D::fillPSetDescription ( edm::ParameterSetDescription iDesc)
static

Definition at line 726 of file PatternRecognitionbyCLUE3D.cc.

References edm::ParameterSetDescription::add().

726  {
727  iDesc.add<int>("algo_verbosity", 0);
728  iDesc.add<std::vector<double>>("criticalDensity", {4, 4, 4})->setComment("in GeV");
729  iDesc.add<std::vector<double>>("criticalSelfDensity", {0.15, 0.15, 0.15} /* roughly 1/(densitySiblingLayers+1) */)
730  ->setComment("Minimum ratio of self_energy/local_density to become a seed.");
731  iDesc.add<std::vector<int>>("densitySiblingLayers", {3, 3, 3})
732  ->setComment(
733  "inclusive, layers to consider while computing local density and searching for nearestHigher higher");
734  iDesc.add<std::vector<double>>("densityEtaPhiDistanceSqr", {0.0008, 0.0008, 0.0008})
735  ->setComment("in eta,phi space, distance to consider for local density");
736  iDesc.add<std::vector<double>>("densityXYDistanceSqr", {3.24, 3.24, 3.24})
737  ->setComment("in cm, distance on the transverse plane to consider for local density");
738  iDesc.add<std::vector<double>>("kernelDensityFactor", {0.2, 0.2, 0.2})
739  ->setComment("Kernel factor to be applied to other LC while computing the local density");
740  iDesc.add<bool>("densityOnSameLayer", false);
741  iDesc.add<bool>("nearestHigherOnSameLayer", false)
742  ->setComment("Allow the nearestHigher to be located on the same layer");
743  iDesc.add<bool>("useAbsoluteProjectiveScale", true)
744  ->setComment("Express all cuts in terms of r/z*z_0{,phi} projective variables");
745  iDesc.add<bool>("useClusterDimensionXY", false)
746  ->setComment(
747  "Boolean. If true use the estimated cluster radius to determine the cluster compatibility while computing "
748  "the local density");
749  iDesc.add<bool>("rescaleDensityByZ", false)
750  ->setComment(
751  "Rescale local density by the extension of the Z 'volume' explored. The transvere dimension is, at present, "
752  "fixed and factored out.");
753  iDesc.add<std::vector<double>>("criticalEtaPhiDistance", {0.025, 0.025, 0.025})
754  ->setComment("Minimal distance in eta,phi space from nearestHigher to become a seed");
755  iDesc.add<std::vector<double>>("criticalXYDistance", {1.8, 1.8, 1.8})
756  ->setComment("Minimal distance in cm on the XY plane from nearestHigher to become a seed");
757  iDesc.add<std::vector<int>>("criticalZDistanceLyr", {5, 5, 5})
758  ->setComment("Minimal distance in layers along the Z axis from nearestHigher to become a seed");
759  iDesc.add<std::vector<double>>("outlierMultiplier", {2, 2, 2})
760  ->setComment("Minimal distance in transverse space from nearestHigher to become an outlier");
761  iDesc.add<std::vector<int>>("minNumLayerCluster", {2, 2, 2})->setComment("Not Inclusive");
762  iDesc.add<bool>("doPidCut", false);
763  iDesc.add<double>("cutHadProb", 0.5);
764  iDesc.add<bool>("computeLocalTime", false);
765  iDesc.add<bool>("usePCACleaning", false)->setComment("Enable PCA cleaning alorithm");
766 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ filter()

template<typename TILES >
void PatternRecognitionbyCLUE3D::filter ( std::vector< Trackster > &  output,
const std::vector< Trackster > &  inTracksters,
const typename PatternRecognitionAlgoBaseT< TILES >::Inputs input,
std::unordered_map< int, std::vector< int >> &  seedToTracksterAssociation 
)
override

Definition at line 359 of file PatternRecognitionbyCLUE3D.cc.

References ticl::Trackster::charged_hadron, ticl::Trackster::neutral_hadron, and submitPVValidationJobs::t.

362  {
363  auto isHAD = [this](const Trackster &t) -> bool {
364  auto const hadProb = t.id_probability(ticl::Trackster::ParticleType::charged_hadron) +
366  return hadProb >= cutHadProb_;
367  };
368 
369  if (doPidCut_) {
370  for (auto const &t : inTracksters) {
371  if (!isHAD(t)) {
372  output.push_back(t);
373  }
374  }
375  } else {
376  output = inTracksters;
377  }
378 }
Definition: output.py:1

◆ findAndAssignTracksters()

template<typename TILES >
int PatternRecognitionbyCLUE3D::findAndAssignTracksters ( const TILES &  tiles,
const std::vector< std::pair< int, int >> &  layerIdx2layerandSoa 
)
private

Definition at line 653 of file PatternRecognitionbyCLUE3D.cc.

References ticl::Advanced, hltPFPuppi_cfi::algoId, mps_fire::i, and electrons_cff::numberOfClusters.

654  {
655  unsigned int nTracksters = 0;
656  std::vector<std::pair<int, int>> localStack;
657  const auto &critical_transverse_distance =
659  // find cluster seeds and outlier
660  for (unsigned int layer = 0; layer < 2 * rhtools_.lastLayer(); layer++) {
661  auto &clustersOnLayer = clusters_[layer];
662  unsigned int numberOfClusters = clustersOnLayer.x.size();
663  for (unsigned int i = 0; i < numberOfClusters; i++) {
664  // initialize clusterIndex
665  clustersOnLayer.clusterIndex[i] = -1;
666  auto algoId = clustersOnLayer.algoId[i];
667  bool isSeed = (clustersOnLayer.delta[i].first > critical_transverse_distance[algoId] ||
668  clustersOnLayer.delta[i].second > criticalZDistanceLyr_[algoId]) &&
669  (clustersOnLayer.rho[i] >= criticalDensity_[algoId]) &&
670  (clustersOnLayer.energy[i] / clustersOnLayer.rho[i] > criticalSelfDensity_[algoId]);
671  if (!clustersOnLayer.isSilicon[i]) {
672  isSeed = (clustersOnLayer.delta[i].first > clustersOnLayer.radius[i] ||
673  clustersOnLayer.delta[i].second > criticalZDistanceLyr_[algoId]) &&
674  (clustersOnLayer.rho[i] >= criticalDensity_[algoId]) &&
675  (clustersOnLayer.energy[i] / clustersOnLayer.rho[i] > criticalSelfDensity_[algoId]);
676  }
677  bool isOutlier =
678  (clustersOnLayer.delta[i].first > outlierMultiplier_[algoId] * critical_transverse_distance[algoId]) &&
679  (clustersOnLayer.rho[i] < criticalDensity_[algoId]);
680  if (isSeed) {
682  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
683  << "Found seed on Layer " << layer << " SOAidx: " << i << " assigned ClusterIdx: " << nTracksters;
684  }
685  clustersOnLayer.clusterIndex[i] = nTracksters++;
686  tracksterSeedAlgoId_.push_back(algoId);
687  clustersOnLayer.isSeed[i] = true;
688  localStack.emplace_back(layer, i);
689  } else if (!isOutlier) {
690  auto [lyrIdx, soaIdx] = clustersOnLayer.nearestHigher[i];
692  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
693  << "Found follower on Layer " << layer << " SOAidx: " << i << " attached to cluster on layer: " << lyrIdx
694  << " SOAidx: " << soaIdx;
695  }
696  if (lyrIdx >= 0)
697  clusters_[lyrIdx].followers[soaIdx].emplace_back(layer, i);
698  } else {
700  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
701  << "Found Outlier on Layer " << layer << " SOAidx: " << i << " with rho: " << clustersOnLayer.rho[i]
702  << " and delta: " << clustersOnLayer.delta[i].first << ", " << clustersOnLayer.delta[i].second;
703  }
704  }
705  }
706  }
707 
708  // Propagate cluster index
709  while (!localStack.empty()) {
710  auto [lyrIdx, soaIdx] = localStack.back();
711  auto &thisSeed = clusters_[lyrIdx].followers[soaIdx];
712  localStack.pop_back();
713 
714  // loop over followers
715  for (auto [follower_lyrIdx, follower_soaIdx] : thisSeed) {
716  // pass id to a follower
717  clusters_[follower_lyrIdx].clusterIndex[follower_soaIdx] = clusters_[lyrIdx].clusterIndex[soaIdx];
718  // push this follower to localStack
719  localStack.emplace_back(follower_lyrIdx, follower_soaIdx);
720  }
721  }
722  return nTracksters;
723 }
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_
const std::vector< double > outlierMultiplier_
const std::vector< double > criticalDensity_
const std::vector< double > criticalXYDistance_
const std::vector< double > criticalEtaPhiDistance_
const std::vector< int > criticalZDistanceLyr_
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:80
const std::vector< double > criticalSelfDensity_

◆ makeTracksters()

template<typename TILES >
void PatternRecognitionbyCLUE3D::makeTracksters ( const typename PatternRecognitionAlgoBaseT< TILES >::Inputs input,
std::vector< Trackster > &  result,
std::unordered_map< int, std::vector< int >> &  seedToTracksterAssociation 
)
override

Definition at line 148 of file PatternRecognitionbyCLUE3D.cc.

References funct::abs(), ticl::Advanced, cms::cuda::assert(), ticl::assignPCAtoTracksters(), hcalRecHitTable_cff::detId, f, relativeConstraints::geom, edm::EventSetup::getData(), reco::CaloCluster::hgcal_em, mps_fire::i, input, createfilelist::int, WZElectronSkims53X_cff::max, point, reset(), mps_fire::result, mathSSE::sqrt(), makeMuonMisalignmentScenario::sum_x, makeMuonMisalignmentScenario::sum_y, submitPVValidationJobs::t, findQualityFiles::v, and x.

151  {
152  // Protect from events with no seeding regions
153  if (input.regions.empty())
154  return;
155 
156  const int eventNumber = input.ev.eventAuxiliary().event();
158  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "New Event";
159  }
160 
161  edm::EventSetup const &es = input.es;
164 
165  // Assume identical Z-positioning between positive and negative sides.
166  // Also, layers inside the HGCAL geometry start from 1.
167  for (unsigned int i = 0; i < rhtools_.lastLayer(); ++i) {
168  layersPosZ_.push_back(rhtools_.getPositionLayer(i + 1).z());
170  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Layer " << i << " located at Z: " << layersPosZ_.back();
171  }
172  }
173 
174  clusters_.clear();
175  tracksterSeedAlgoId_.clear();
176 
177  clusters_.resize(2 * rhtools_.lastLayer(false));
178  std::vector<std::pair<int, int>> layerIdx2layerandSoa; //used everywhere also to propagate cluster masking
179 
180  layerIdx2layerandSoa.reserve(input.layerClusters.size());
181  unsigned int layerIdx = 0;
182  for (auto const &lc : input.layerClusters) {
183  if (input.mask[layerIdx] == 0.) {
185  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Skipping masked cluster: " << layerIdx;
186  }
187  layerIdx2layerandSoa.emplace_back(-1, -1);
188  layerIdx++;
189  continue;
190  }
191  const auto firstHitDetId = lc.hitsAndFractions()[0].first;
192  int layer = rhtools_.getLayerWithOffset(firstHitDetId) - 1 +
193  rhtools_.lastLayer(false) * ((rhtools_.zside(firstHitDetId) + 1) >> 1);
194  assert(layer >= 0);
195  auto detId = lc.hitsAndFractions()[0].first;
196  int layerClusterIndexInLayer = clusters_[layer].x.size();
197  layerIdx2layerandSoa.emplace_back(layer, layerClusterIndexInLayer);
198  float sum_x = 0.;
199  float sum_y = 0.;
200  float sum_sqr_x = 0.;
201  float sum_sqr_y = 0.;
202  float ref_x = lc.x();
203  float ref_y = lc.y();
204  float invClsize = 1. / lc.hitsAndFractions().size();
205  for (auto const &hitsAndFractions : lc.hitsAndFractions()) {
206  auto const &point = rhtools_.getPosition(hitsAndFractions.first);
207  sum_x += point.x() - ref_x;
208  sum_sqr_x += (point.x() - ref_x) * (point.x() - ref_x);
209  sum_y += point.y() - ref_y;
210  sum_sqr_y += (point.y() - ref_y) * (point.y() - ref_y);
211  }
212  // The variance of X for X uniform in circle of radius R is R^2/4,
213  // therefore we multiply the sqrt(var) by 2 to have a rough estimate of the
214  // radius. On the other hand, while averaging the x and y radius, we would
215  // end up dividing by 2. Hence we omit the value here and in the average
216  // below, too.
217  float radius_x = sqrt((sum_sqr_x - (sum_x * sum_x) * invClsize) * invClsize);
218  float radius_y = sqrt((sum_sqr_y - (sum_y * sum_y) * invClsize) * invClsize);
220  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
221  << "cluster rx: " << std::setw(5) << radius_x << ", ry: " << std::setw(5) << radius_y
222  << ", r: " << std::setw(5) << (radius_x + radius_y) << ", cells: " << std::setw(4)
223  << lc.hitsAndFractions().size();
224  }
225 
226  // The case of single cell layer clusters has to be handled differently.
227 
228  if (invClsize == 1.) {
229  // Silicon case
230  if (rhtools_.isSilicon(detId)) {
231  radius_x = radius_y = rhtools_.getRadiusToSide(detId);
233  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Single cell cluster in silicon, rx: " << std::setw(5)
234  << radius_x << ", ry: " << std::setw(5) << radius_y;
235  }
236  } else {
237  auto const &point = rhtools_.getPosition(detId);
238  auto const &eta_phi_window = rhtools_.getScintDEtaDPhi(detId);
239  radius_x = radius_y = point.perp() * eta_phi_window.second;
241  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
242  << "Single cell cluster in scintillator. rx: " << std::setw(5) << radius_x << ", ry: " << std::setw(5)
243  << radius_y << ", eta-span: " << std::setw(5) << eta_phi_window.first << ", phi-span: " << std::setw(5)
244  << eta_phi_window.second;
245  }
246  }
247  }
248 
249  // Maybe check if these vectors can be reserved beforehands
250  clusters_[layer].x.emplace_back(lc.x());
251  clusters_[layer].y.emplace_back(lc.y());
252  clusters_[layer].z.emplace_back(lc.z());
253  clusters_[layer].r_over_absz.emplace_back(sqrt(lc.x() * lc.x() + lc.y() * lc.y()) / std::abs(lc.z()));
254  clusters_[layer].radius.emplace_back(radius_x + radius_y);
255  clusters_[layer].eta.emplace_back(lc.eta());
256  clusters_[layer].phi.emplace_back(lc.phi());
257  clusters_[layer].cells.push_back(lc.hitsAndFractions().size());
258  clusters_[layer].algoId.push_back(lc.algo() - reco::CaloCluster::hgcal_em);
259  clusters_[layer].isSilicon.push_back(rhtools_.isSilicon(detId));
260  clusters_[layer].energy.emplace_back(lc.energy());
261  clusters_[layer].isSeed.push_back(false);
262  clusters_[layer].clusterIndex.emplace_back(-1);
263  clusters_[layer].layerClusterOriginalIdx.emplace_back(layerIdx++);
264  clusters_[layer].nearestHigher.emplace_back(-1, -1);
265  clusters_[layer].rho.emplace_back(0.f);
266  clusters_[layer].z_extension.emplace_back(0.f);
267  clusters_[layer].delta.emplace_back(
269  }
270  for (unsigned int layer = 0; layer < clusters_.size(); layer++) {
271  clusters_[layer].followers.resize(clusters_[layer].x.size());
272  }
273 
274  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
275  int maxLayer = 2 * lastLayerPerSide - 1;
276  std::vector<int> numberOfClustersPerLayer(maxLayer, 0);
277  for (int i = 0; i <= maxLayer; i++) {
278  calculateLocalDensity(input.tiles, i, layerIdx2layerandSoa);
279  }
280  for (int i = 0; i <= maxLayer; i++) {
281  calculateDistanceToHigher(input.tiles, i, layerIdx2layerandSoa);
282  }
283 
284  auto nTracksters = findAndAssignTracksters(input.tiles, layerIdx2layerandSoa);
286  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Reconstructed " << nTracksters << " tracksters" << std::endl;
287  dumpClusters(input.tiles, layerIdx2layerandSoa, eventNumber);
288  }
289 
290  // Build Trackster
291  result.resize(nTracksters);
292 
293  for (unsigned int layer = 0; layer < clusters_.size(); ++layer) {
294  const auto &thisLayer = clusters_[layer];
296  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Examining Layer: " << layer;
297  }
298  for (unsigned int lc = 0; lc < thisLayer.x.size(); ++lc) {
300  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Trackster " << thisLayer.clusterIndex[lc];
301  }
302  if (thisLayer.clusterIndex[lc] >= 0) {
304  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << " adding lcIdx: " << thisLayer.layerClusterOriginalIdx[lc];
305  }
306  result[thisLayer.clusterIndex[lc]].vertices().push_back(thisLayer.layerClusterOriginalIdx[lc]);
307  result[thisLayer.clusterIndex[lc]].vertex_multiplicity().push_back(1);
308  // loop over followers
309  for (auto [follower_lyrIdx, follower_soaIdx] : thisLayer.followers[lc]) {
310  std::array<unsigned int, 2> edge = {
311  {(unsigned int)thisLayer.layerClusterOriginalIdx[lc],
312  (unsigned int)clusters_[follower_lyrIdx].layerClusterOriginalIdx[follower_soaIdx]}};
313  result[thisLayer.clusterIndex[lc]].edges().push_back(edge);
314  }
315  }
316  }
317  }
318  size_t tracksterIndex = 0;
319  result.erase(std::remove_if(std::begin(result),
320  std::end(result),
321  [&](auto const &v) {
322  return static_cast<int>(v.vertices().size()) <
323  minNumLayerCluster_[tracksterSeedAlgoId_[tracksterIndex++]];
324  }),
325  result.end());
326 
327  result.shrink_to_fit();
328 
330  input.layerClusters,
331  input.layerClustersTime,
333  rhtools_,
335  true, // energy weighting
337 
339  for (auto const &t : result) {
340  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Barycenter: " << t.barycenter();
341  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "LCs: " << t.vertices().size();
342  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Energy: " << t.raw_energy();
343  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Regressed: " << t.regressed_energy();
344  }
345  }
346 
347  // Dump Tracksters information
349  dumpTracksters(layerIdx2layerandSoa, eventNumber, result);
350  }
351 
352  // Reset internal clusters_ structure of array for next event
353  reset();
355  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << std::endl;
356  }
357 }
Log< level::Info, true > LogVerbatim
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
T z() const
Definition: PV3DBase.h:61
std::float_t getRadiusToSide(const DetId &) const
Definition: RecHitTools.cc:246
void assignPCAtoTracksters(std::vector< Trackster > &tracksters, const std::vector< reco::CaloCluster > &layerClusters, const edm::ValueMap< std::pair< float, float >> &layerClustersTime, double z_limit_em, hgcal::RecHitTools const &rhTools, bool computeLocalTime=false, bool energyWeight=true, bool clean=false, int minLayer=10, int maxLayer=10)
std::vector< ClustersOnLayer > clusters_
std::pair< float, float > getScintDEtaDPhi(const DetId &) const
Definition: RecHitTools.cc:236
assert(be >=bs)
int zside(const DetId &id) const
Definition: RecHitTools.cc:174
void calculateDistanceToHigher(const TILES &, const int layerId, const std::vector< std::pair< int, int >> &)
static std::string const input
Definition: EdmProvDump.cc:50
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:140
T sqrt(T t)
Definition: SSEVec.h:23
bool isSilicon(const DetId &) const
Definition: RecHitTools.cc:478
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
int findAndAssignTracksters(const TILES &, const std::vector< std::pair< int, int >> &)
void dumpClusters(const TILES &tiles, const std::vector< std::pair< int, int >> &layerIdx2layerandSoa, const int) const
GlobalPoint getPositionLayer(int layer, bool nose=false) const
Definition: RecHitTools.cc:152
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:79
void calculateLocalDensity(const TILES &, const int layerId, const std::vector< std::pair< int, int >> &)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
unsigned int lastLayerEE(bool nose=false) const
Definition: RecHitTools.h:76
void dumpTracksters(const std::vector< std::pair< int, int >> &layerIdx2layerandSoa, const int, const std::vector< Trackster > &) const
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:80
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381

◆ reset()

template<typename TILES >
void ticl::PatternRecognitionbyCLUE3D< TILES >::reset ( void  )
inlineprivate

Definition at line 101 of file PatternRecognitionbyCLUE3D.h.

References DummyCfis::c, and ticl::PatternRecognitionbyCLUE3D< TILES >::clusters_.

101  {
102  for (auto& c : clusters_) {
103  c.clear();
104  c.shrink_to_fit();
105  }
106  }
std::vector< ClustersOnLayer > clusters_

Member Data Documentation

◆ caloGeomToken_

template<typename TILES >
edm::ESGetToken<CaloGeometry, CaloGeometryRecord> ticl::PatternRecognitionbyCLUE3D< TILES >::caloGeomToken_
private

Definition at line 122 of file PatternRecognitionbyCLUE3D.h.

◆ clusters_

template<typename TILES >
std::vector<ClustersOnLayer> ticl::PatternRecognitionbyCLUE3D< TILES >::clusters_
private

◆ computeLocalTime_

template<typename TILES >
const bool ticl::PatternRecognitionbyCLUE3D< TILES >::computeLocalTime_
private

Definition at line 142 of file PatternRecognitionbyCLUE3D.h.

◆ criticalDensity_

template<typename TILES >
const std::vector<double> ticl::PatternRecognitionbyCLUE3D< TILES >::criticalDensity_
private

Definition at line 123 of file PatternRecognitionbyCLUE3D.h.

◆ criticalEtaPhiDistance_

template<typename TILES >
const std::vector<double> ticl::PatternRecognitionbyCLUE3D< TILES >::criticalEtaPhiDistance_
private

Definition at line 134 of file PatternRecognitionbyCLUE3D.h.

◆ criticalSelfDensity_

template<typename TILES >
const std::vector<double> ticl::PatternRecognitionbyCLUE3D< TILES >::criticalSelfDensity_
private

Definition at line 124 of file PatternRecognitionbyCLUE3D.h.

◆ criticalXYDistance_

template<typename TILES >
const std::vector<double> ticl::PatternRecognitionbyCLUE3D< TILES >::criticalXYDistance_
private

Definition at line 135 of file PatternRecognitionbyCLUE3D.h.

◆ criticalZDistanceLyr_

template<typename TILES >
const std::vector<int> ticl::PatternRecognitionbyCLUE3D< TILES >::criticalZDistanceLyr_
private

Definition at line 136 of file PatternRecognitionbyCLUE3D.h.

◆ cutHadProb_

template<typename TILES >
const float ticl::PatternRecognitionbyCLUE3D< TILES >::cutHadProb_
private

Definition at line 140 of file PatternRecognitionbyCLUE3D.h.

◆ densityEtaPhiDistanceSqr_

template<typename TILES >
const std::vector<double> ticl::PatternRecognitionbyCLUE3D< TILES >::densityEtaPhiDistanceSqr_
private

Definition at line 126 of file PatternRecognitionbyCLUE3D.h.

◆ densityOnSameLayer_

template<typename TILES >
const bool ticl::PatternRecognitionbyCLUE3D< TILES >::densityOnSameLayer_
private

Definition at line 129 of file PatternRecognitionbyCLUE3D.h.

◆ densitySiblingLayers_

template<typename TILES >
const std::vector<int> ticl::PatternRecognitionbyCLUE3D< TILES >::densitySiblingLayers_
private

Definition at line 125 of file PatternRecognitionbyCLUE3D.h.

◆ densityXYDistanceSqr_

template<typename TILES >
const std::vector<double> ticl::PatternRecognitionbyCLUE3D< TILES >::densityXYDistanceSqr_
private

Definition at line 127 of file PatternRecognitionbyCLUE3D.h.

◆ doPidCut_

template<typename TILES >
const bool ticl::PatternRecognitionbyCLUE3D< TILES >::doPidCut_
private

Definition at line 139 of file PatternRecognitionbyCLUE3D.h.

◆ filter_on_categories_

template<typename TILES >
const std::vector<int> ticl::PatternRecognitionbyCLUE3D< TILES >::filter_on_categories_
private

Definition at line 141 of file PatternRecognitionbyCLUE3D.h.

◆ kernelDensityFactor_

template<typename TILES >
const std::vector<double> ticl::PatternRecognitionbyCLUE3D< TILES >::kernelDensityFactor_
private

Definition at line 128 of file PatternRecognitionbyCLUE3D.h.

◆ layersPosZ_

template<typename TILES >
std::vector<float> ticl::PatternRecognitionbyCLUE3D< TILES >::layersPosZ_
private

Definition at line 119 of file PatternRecognitionbyCLUE3D.h.

◆ minNumLayerCluster_

template<typename TILES >
const std::vector<int> ticl::PatternRecognitionbyCLUE3D< TILES >::minNumLayerCluster_
private

Definition at line 138 of file PatternRecognitionbyCLUE3D.h.

◆ nearestHigherOnSameLayer_

template<typename TILES >
const bool ticl::PatternRecognitionbyCLUE3D< TILES >::nearestHigherOnSameLayer_
private

Definition at line 130 of file PatternRecognitionbyCLUE3D.h.

◆ outlierMultiplier_

template<typename TILES >
const std::vector<double> ticl::PatternRecognitionbyCLUE3D< TILES >::outlierMultiplier_
private

Definition at line 137 of file PatternRecognitionbyCLUE3D.h.

◆ rescaleDensityByZ_

template<typename TILES >
const bool ticl::PatternRecognitionbyCLUE3D< TILES >::rescaleDensityByZ_
private

Definition at line 133 of file PatternRecognitionbyCLUE3D.h.

◆ rhtools_

template<typename TILES >
hgcal::RecHitTools ticl::PatternRecognitionbyCLUE3D< TILES >::rhtools_
private

Definition at line 145 of file PatternRecognitionbyCLUE3D.h.

◆ tracksterSeedAlgoId_

template<typename TILES >
std::vector<int> ticl::PatternRecognitionbyCLUE3D< TILES >::tracksterSeedAlgoId_
private

Definition at line 120 of file PatternRecognitionbyCLUE3D.h.

◆ useAbsoluteProjectiveScale_

template<typename TILES >
const bool ticl::PatternRecognitionbyCLUE3D< TILES >::useAbsoluteProjectiveScale_
private

Definition at line 131 of file PatternRecognitionbyCLUE3D.h.

◆ useClusterDimensionXY_

template<typename TILES >
const bool ticl::PatternRecognitionbyCLUE3D< TILES >::useClusterDimensionXY_
private

Definition at line 132 of file PatternRecognitionbyCLUE3D.h.

◆ usePCACleaning_

template<typename TILES >
const bool ticl::PatternRecognitionbyCLUE3D< TILES >::usePCACleaning_
private

Definition at line 143 of file PatternRecognitionbyCLUE3D.h.