CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static 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 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 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::string eidInputName_
 
const float eidMinClusterEnergy_
 
const int eidNClusters_
 
const int eidNLayers_
 
const std::string eidOutputNameEnergy_
 
const std::string eidOutputNameId_
 
tensorflow::Session * eidSession_
 
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_
 

Static Private Attributes

static const int eidNFeatures_ = 3
 

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  eidInputName_(conf.getParameter<std::string>("eid_input_name")),
43  eidOutputNameEnergy_(conf.getParameter<std::string>("eid_output_name_energy")),
44  eidOutputNameId_(conf.getParameter<std::string>("eid_output_name_id")),
45  eidMinClusterEnergy_(conf.getParameter<double>("eid_min_cluster_energy")),
46  eidNLayers_(conf.getParameter<int>("eid_n_layers")),
47  eidNClusters_(conf.getParameter<int>("eid_n_clusters")),
48  computeLocalTime_(conf.getParameter<bool>("computeLocalTime")){};
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 680 of file PatternRecognitionbyCLUE3D.cc.

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

681  {
684  auto &clustersOnLayer = clusters_[layerId];
685  unsigned int numberOfClusters = clustersOnLayer.x.size();
686 
687  auto distanceSqr = [](float r0, float r1, float phi0, float phi1) -> float {
688  auto delta_phi = reco::deltaPhi(phi0, phi1);
689  return (r0 - r1) * (r0 - r1) + r1 * r1 * delta_phi * delta_phi;
690  };
691 
692  for (unsigned int i = 0; i < numberOfClusters; i++) {
694  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
695  << "Starting searching nearestHigher on " << layerId << " with rho: " << clustersOnLayer.rho[i]
696  << " at eta, phi: " << tiles[layerId].etaBin(clustersOnLayer.eta[i]) << ", "
697  << tiles[layerId].phiBin(clustersOnLayer.phi[i]);
698  }
699  // We need to partition the two sides of the HGCAL detector
700  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
701  int minLayer = 0;
702  int maxLayer = 2 * lastLayerPerSide - 1;
703  auto algoId = clustersOnLayer.algoId[i];
704  if (layerId < lastLayerPerSide) {
706  maxLayer = std::min(layerId + densitySiblingLayers_[algoId], lastLayerPerSide - 1);
707  } else {
708  minLayer = std::max(layerId - densitySiblingLayers_[algoId], lastLayerPerSide + 1);
709  maxLayer = std::min(layerId + densitySiblingLayers_[algoId], maxLayer);
710  }
712  float i_delta = maxDelta;
713  std::pair<int, int> i_nearestHigher(-1, -1);
714  std::pair<float, int> nearest_distances(maxDelta, std::numeric_limits<int>::max());
715  for (int currentLayer = minLayer; currentLayer <= maxLayer; currentLayer++) {
716  if (!nearestHigherOnSameLayer_ && (layerId == currentLayer))
717  continue;
718  const auto &tileOnLayer = tiles[currentLayer];
719  int etaWindow = 1;
720  int phiWindow = 1;
721  int etaBinMin = std::max(tileOnLayer.etaBin(clustersOnLayer.eta[i]) - etaWindow, 0);
722  int etaBinMax = std::min(tileOnLayer.etaBin(clustersOnLayer.eta[i]) + etaWindow, nEtaBin - 1);
723  int phiBinMin = tileOnLayer.phiBin(clustersOnLayer.phi[i]) - phiWindow;
724  int phiBinMax = tileOnLayer.phiBin(clustersOnLayer.phi[i]) + phiWindow;
725  for (int ieta = etaBinMin; ieta <= etaBinMax; ++ieta) {
726  auto offset = ieta * nPhiBin;
727  for (int iphi_it = phiBinMin; iphi_it <= phiBinMax; ++iphi_it) {
728  int iphi = ((iphi_it % nPhiBin + nPhiBin) % nPhiBin);
730  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
731  << "Searching nearestHigher on " << currentLayer << " eta, phi: " << ieta << ", " << iphi_it << " "
732  << iphi << " " << offset << " " << (offset + iphi);
733  }
734  for (auto otherClusterIdx : tileOnLayer[offset + iphi]) {
735  auto const &layerandSoa = layerIdx2layerandSoa[otherClusterIdx];
736  // Skip masked layer clusters
737  if ((layerandSoa.first == -1) && (layerandSoa.second == -1))
738  continue;
739  auto const &clustersOnOtherLayer = clusters_[layerandSoa.first];
740  auto dist = maxDelta;
741  auto dist_transverse = maxDelta;
742  int dist_layers = std::abs(layerandSoa.first - layerId);
744  dist_transverse = distanceSqr(clustersOnLayer.r_over_absz[i] * clustersOnLayer.z[i],
745  clustersOnOtherLayer.r_over_absz[layerandSoa.second] * clustersOnLayer.z[i],
746  clustersOnLayer.phi[i],
747  clustersOnOtherLayer.phi[layerandSoa.second]);
748  // Add Z-scale to the final distance
749  dist = dist_transverse;
750  } else {
751  dist = reco::deltaR2(clustersOnLayer.eta[i],
752  clustersOnLayer.phi[i],
753  clustersOnOtherLayer.eta[layerandSoa.second],
754  clustersOnOtherLayer.phi[layerandSoa.second]);
755  dist_transverse = dist;
756  }
757  bool foundHigher = (clustersOnOtherLayer.rho[layerandSoa.second] > clustersOnLayer.rho[i]) ||
758  (clustersOnOtherLayer.rho[layerandSoa.second] == clustersOnLayer.rho[i] &&
759  clustersOnOtherLayer.layerClusterOriginalIdx[layerandSoa.second] >
760  clustersOnLayer.layerClusterOriginalIdx[i]);
762  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
763  << "Searching nearestHigher on " << currentLayer
764  << " with rho: " << clustersOnOtherLayer.rho[layerandSoa.second]
765  << " on layerIdxInSOA: " << layerandSoa.first << ", " << layerandSoa.second
766  << " with distance: " << sqrt(dist) << " foundHigher: " << foundHigher;
767  }
768  if (foundHigher && dist <= i_delta) {
769  // update i_delta
770  i_delta = dist;
771  nearest_distances = std::make_pair(sqrt(dist_transverse), dist_layers);
772  // update i_nearestHigher
773  i_nearestHigher = layerandSoa;
774  }
775  } // End of loop on clusters
776  } // End of loop on phi bins
777  } // End of loop on eta bins
778  } // End of loop on layers
779 
780  clustersOnLayer.delta[i] = nearest_distances;
781  clustersOnLayer.nearestHigher[i] = i_nearestHigher;
782  } // End of loop on clusters
783 }
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 514 of file PatternRecognitionbyCLUE3D.cc.

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

515  {
518  auto &clustersOnLayer = clusters_[layerId];
519  unsigned int numberOfClusters = clustersOnLayer.x.size();
520 
521  auto isReachable = [](float r0, float r1, float phi0, float phi1, float delta_sqr) -> bool {
522  auto delta_phi = reco::deltaPhi(phi0, phi1);
523  return (r0 - r1) * (r0 - r1) + r1 * r1 * delta_phi * delta_phi < delta_sqr;
524  };
525  auto distance_debug = [&](float x1, float x2, float y1, float y2) -> float {
526  return sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2));
527  };
528 
529  for (unsigned int i = 0; i < numberOfClusters; i++) {
530  auto algoId = clustersOnLayer.algoId[i];
531  // We need to partition the two sides of the HGCAL detector
532  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
533  int minLayer = 0;
534  int maxLayer = 2 * lastLayerPerSide - 1;
535  if (layerId < lastLayerPerSide) {
537  maxLayer = std::min(layerId + densitySiblingLayers_[algoId], lastLayerPerSide - 1);
538  } else {
539  minLayer = std::max(layerId - densitySiblingLayers_[algoId], lastLayerPerSide);
540  maxLayer = std::min(layerId + densitySiblingLayers_[algoId], maxLayer);
541  }
542  float deltaLayersZ = std::abs(layersPosZ_[maxLayer % lastLayerPerSide] - layersPosZ_[minLayer % lastLayerPerSide]);
543 
544  for (int currentLayer = minLayer; currentLayer <= maxLayer; currentLayer++) {
546  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "RefLayer: " << layerId << " SoaIDX: " << i;
547  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "NextLayer: " << currentLayer;
548  }
549  const auto &tileOnLayer = tiles[currentLayer];
550  bool onSameLayer = (currentLayer == layerId);
552  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "onSameLayer: " << onSameLayer;
553  }
554  const int etaWindow = 2;
555  const int phiWindow = 2;
556  int etaBinMin = std::max(tileOnLayer.etaBin(clustersOnLayer.eta[i]) - etaWindow, 0);
557  int etaBinMax = std::min(tileOnLayer.etaBin(clustersOnLayer.eta[i]) + etaWindow, nEtaBin - 1);
558  int phiBinMin = tileOnLayer.phiBin(clustersOnLayer.phi[i]) - phiWindow;
559  int phiBinMax = tileOnLayer.phiBin(clustersOnLayer.phi[i]) + phiWindow;
561  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "eta: " << clustersOnLayer.eta[i];
562  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "phi: " << clustersOnLayer.phi[i];
563  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "etaBinMin: " << etaBinMin << ", etaBinMax: " << etaBinMax;
564  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "phiBinMin: " << phiBinMin << ", phiBinMax: " << phiBinMax;
565  }
566  for (int ieta = etaBinMin; ieta <= etaBinMax; ++ieta) {
567  auto offset = ieta * nPhiBin;
569  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "offset: " << offset;
570  }
571  for (int iphi_it = phiBinMin; iphi_it <= phiBinMax; ++iphi_it) {
572  int iphi = ((iphi_it % nPhiBin + nPhiBin) % nPhiBin);
574  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "iphi: " << iphi;
575  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
576  << "Entries in tileBin: " << tileOnLayer[offset + iphi].size();
577  }
578  for (auto otherClusterIdx : tileOnLayer[offset + iphi]) {
579  auto const &layerandSoa = layerIdx2layerandSoa[otherClusterIdx];
580  // Skip masked layer clusters
581  if ((layerandSoa.first == -1) && (layerandSoa.second == -1)) {
583  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Skipping masked layerIdx " << otherClusterIdx;
584  }
585  continue;
586  }
587  auto const &clustersLayer = clusters_[layerandSoa.first];
589  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
590  << "OtherLayer: " << layerandSoa.first << " SoaIDX: " << layerandSoa.second;
591  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "OtherEta: " << clustersLayer.eta[layerandSoa.second];
592  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "OtherPhi: " << clustersLayer.phi[layerandSoa.second];
593  }
594 
595  bool onSameCluster = clustersOnLayer.layerClusterOriginalIdx[i] == otherClusterIdx;
596  if (onSameLayer && !densityOnSameLayer_ && !onSameCluster) {
598  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
599  << "Skipping different cluster " << otherClusterIdx << "in the same layer " << currentLayer;
600  }
601  continue;
602  }
603 
604  bool reachable = false;
607  reachable = isReachable(clustersOnLayer.r_over_absz[i] * clustersOnLayer.z[i],
608  clustersLayer.r_over_absz[layerandSoa.second] * clustersOnLayer.z[i],
609  clustersOnLayer.phi[i],
610  clustersLayer.phi[layerandSoa.second],
611  clustersOnLayer.radius[i] * clustersOnLayer.radius[i]);
612  } else {
613  // Still differentiate between silicon and Scintillator.
614  // Scintillator has yet to be studied further.
615  if (clustersOnLayer.isSilicon[i]) {
616  reachable = isReachable(clustersOnLayer.r_over_absz[i] * clustersOnLayer.z[i],
617  clustersLayer.r_over_absz[layerandSoa.second] * clustersOnLayer.z[i],
618  clustersOnLayer.phi[i],
619  clustersLayer.phi[layerandSoa.second],
621  } else {
622  reachable = isReachable(clustersOnLayer.r_over_absz[i] * clustersOnLayer.z[i],
623  clustersLayer.r_over_absz[layerandSoa.second] * clustersOnLayer.z[i],
624  clustersOnLayer.phi[i],
625  clustersLayer.phi[layerandSoa.second],
626  clustersOnLayer.radius[i] * clustersOnLayer.radius[i]);
627  }
628  }
629  } else {
630  reachable = (reco::deltaR2(clustersOnLayer.eta[i],
631  clustersOnLayer.phi[i],
632  clustersLayer.eta[layerandSoa.second],
633  clustersLayer.phi[layerandSoa.second]) < densityEtaPhiDistanceSqr_[algoId]);
634  }
636  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Distance[eta,phi]: "
637  << reco::deltaR2(clustersOnLayer.eta[i],
638  clustersOnLayer.phi[i],
639  clustersLayer.eta[layerandSoa.second],
640  clustersLayer.phi[layerandSoa.second]);
641  auto dist = distance_debug(
642  clustersOnLayer.r_over_absz[i],
643  clustersLayer.r_over_absz[layerandSoa.second],
644  clustersOnLayer.r_over_absz[i] * std::abs(clustersOnLayer.phi[i]),
645  clustersLayer.r_over_absz[layerandSoa.second] * std::abs(clustersLayer.phi[layerandSoa.second]));
646  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Distance[cm]: " << (dist * clustersOnLayer.z[i]);
647  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
648  << "Energy Other: " << clustersLayer.energy[layerandSoa.second];
649  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Cluster radius: " << clustersOnLayer.radius[i];
650  }
651  if (reachable) {
652  float factor_same_layer_different_cluster = (onSameLayer && !densityOnSameLayer_) ? 0.f : 1.f;
653  auto energyToAdd = (clustersOnLayer.layerClusterOriginalIdx[i] == otherClusterIdx
654  ? 1.f
655  : kernelDensityFactor_[algoId] * factor_same_layer_different_cluster) *
656  clustersLayer.energy[layerandSoa.second];
657  clustersOnLayer.rho[i] += energyToAdd;
658  clustersOnLayer.z_extension[i] = deltaLayersZ;
660  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
661  << "Adding " << energyToAdd << " partial " << clustersOnLayer.rho[i];
662  }
663  }
664  } // end of loop on possible compatible clusters
665  } // end of loop over phi-bin region
666  } // end of loop over eta-bin region
667  } // end of loop on the sibling layers
668  if (rescaleDensityByZ_) {
670  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
671  << "Rescaling original density: " << clustersOnLayer.rho[i] << " by Z: " << deltaLayersZ
672  << " to final density/cm: " << clustersOnLayer.rho[i] / deltaLayersZ;
673  }
674  clustersOnLayer.rho[i] /= deltaLayersZ;
675  }
676  } // end of loop over clusters on this layer
677 }
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 108 of file PatternRecognitionbyCLUE3D.cc.

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

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

◆ dumpTiles()

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

Definition at line 51 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 hltrates_dqm_sourceclient-live_cfg::offset.

51  {
54  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
55  int maxLayer = 2 * lastLayerPerSide - 1;
56  for (int layer = 0; layer <= maxLayer; layer++) {
57  for (int ieta = 0; ieta < nEtaBin; ieta++) {
58  auto offset = ieta * nPhiBin;
59  for (int phi = 0; phi < nPhiBin; phi++) {
60  int iphi = ((phi % nPhiBin + nPhiBin) % nPhiBin);
61  if (!tiles[layer][offset + iphi].empty()) {
63  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Layer: " << layer << " ieta: " << ieta << " phi: " << phi
64  << " " << tiles[layer][offset + iphi].size();
65  }
66  }
67  }
68  }
69  }
70 }
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 73 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, and findQualityFiles::v.

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

◆ energyRegressionAndID()

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

Definition at line 374 of file PatternRecognitionbyCLUE3D.cc.

References a, funct::abs(), b, HLT_FULL_cff::batchSize, data, hcalRecHitTable_cff::energy, reco::CaloCluster::energy(), reco::CaloCluster::eta(), f, lowptgsfeleseed::features(), nano_mu_digi_cff::float, reco::CaloCluster::hitsAndFractions(), mps_fire::i, input, crabTemplate::inputList, createfilelist::int, dqmiolumiharvest::j, dqmdumpme::k, MainPageGenerator::l, hltEgammaHGCALIDVarsL1Seeded_cfi::layerClusters, jetsAK4_CHS_cff::outputNames, PatBasicFWLiteJetAnalyzer_Selector_cfg::outputs, reco::CaloCluster::phi(), tensorflow::run(), l1trig_cff::shape, jetUpdater_cfi::sort, bphysicsOniaDQM_cfi::vertex, ticl::Trackster::vertex_multiplicity(), AlignmentTracksFromVertexSelector_cfi::vertices, and ticl::Trackster::vertices().

376  {
377  // Energy regression and particle identification strategy:
378  //
379  // 1. Set default values for regressed energy and particle id for each trackster.
380  // 2. Store indices of tracksters whose total sum of cluster energies is above the
381  // eidMinClusterEnergy_ (GeV) threshold. Inference is not applied for soft tracksters.
382  // 3. When no trackster passes the selection, return.
383  // 4. Create input and output tensors. The batch dimension is determined by the number of
384  // selected tracksters.
385  // 5. Fill input tensors with layer cluster features. Per layer, clusters are ordered descending
386  // by energy. Given that tensor data is contiguous in memory, we can use pointer arithmetic to
387  // fill values, even with batching.
388  // 6. Zero-fill features for empty clusters in each layer.
389  // 7. Batched inference.
390  // 8. Assign the regressed energy and id probabilities to each trackster.
391  //
392  // Indices used throughout this method:
393  // i -> batch element / trackster
394  // j -> layer
395  // k -> cluster
396  // l -> feature
397 
398  // set default values per trackster, determine if the cluster energy threshold is passed,
399  // and store indices of hard tracksters
400  std::vector<int> tracksterIndices;
401  for (int i = 0; i < static_cast<int>(tracksters.size()); i++) {
402  // calculate the cluster energy sum (2)
403  // note: after the loop, sumClusterEnergy might be just above the threshold which is enough to
404  // decide whether to run inference for the trackster or not
405  float sumClusterEnergy = 0.;
406  for (const unsigned int &vertex : tracksters[i].vertices()) {
407  sumClusterEnergy += static_cast<float>(layerClusters[vertex].energy());
408  // there might be many clusters, so try to stop early
409  if (sumClusterEnergy >= eidMinClusterEnergy_) {
410  // set default values (1)
411  tracksters[i].setRegressedEnergy(0.f);
412  tracksters[i].zeroProbabilities();
413  tracksterIndices.push_back(i);
414  break;
415  }
416  }
417  }
418 
419  // do nothing when no trackster passes the selection (3)
420  int batchSize = static_cast<int>(tracksterIndices.size());
421  if (batchSize == 0) {
422  return;
423  }
424 
425  // create input and output tensors (4)
426  tensorflow::TensorShape shape({batchSize, eidNLayers_, eidNClusters_, eidNFeatures_});
427  tensorflow::Tensor input(tensorflow::DT_FLOAT, shape);
429 
430  std::vector<tensorflow::Tensor> outputs;
431  std::vector<std::string> outputNames;
432  if (!eidOutputNameEnergy_.empty()) {
434  }
435  if (!eidOutputNameId_.empty()) {
436  outputNames.push_back(eidOutputNameId_);
437  }
438 
439  // fill input tensor (5)
440  for (int i = 0; i < batchSize; i++) {
441  const Trackster &trackster = tracksters[tracksterIndices[i]];
442 
443  // per layer, we only consider the first eidNClusters_ clusters in terms of energy, so in order
444  // to avoid creating large / nested structures to do the sorting for an unknown number of total
445  // clusters, create a sorted list of layer cluster indices to keep track of the filled clusters
446  std::vector<int> clusterIndices(trackster.vertices().size());
447  for (int k = 0; k < (int)trackster.vertices().size(); k++) {
448  clusterIndices[k] = k;
449  }
450  sort(clusterIndices.begin(), clusterIndices.end(), [&layerClusters, &trackster](const int &a, const int &b) {
451  return layerClusters[trackster.vertices(a)].energy() > layerClusters[trackster.vertices(b)].energy();
452  });
453 
454  // keep track of the number of seen clusters per layer
455  std::vector<int> seenClusters(eidNLayers_);
456 
457  // loop through clusters by descending energy
458  for (const int &k : clusterIndices) {
459  // get features per layer and cluster and store the values directly in the input tensor
460  const reco::CaloCluster &cluster = layerClusters[trackster.vertices(k)];
461  int j = rhtools_.getLayerWithOffset(cluster.hitsAndFractions()[0].first) - 1;
462  if (j < eidNLayers_ && seenClusters[j] < eidNClusters_) {
463  // get the pointer to the first feature value for the current batch, layer and cluster
464  float *features = &input.tensor<float, 4>()(i, j, seenClusters[j], 0);
465 
466  // fill features
467  *(features++) = float(cluster.energy() / float(trackster.vertex_multiplicity(k)));
468  *(features++) = float(std::abs(cluster.eta()));
469  *(features) = float(cluster.phi());
470 
471  // increment seen clusters
472  seenClusters[j]++;
473  }
474  }
475 
476  // zero-fill features of empty clusters in each layer (6)
477  for (int j = 0; j < eidNLayers_; j++) {
478  for (int k = seenClusters[j]; k < eidNClusters_; k++) {
479  float *features = &input.tensor<float, 4>()(i, j, k, 0);
480  for (int l = 0; l < eidNFeatures_; l++) {
481  *(features++) = 0.f;
482  }
483  }
484  }
485  }
486 
487  // run the inference (7)
489 
490  // store regressed energy per trackster (8)
491  if (!eidOutputNameEnergy_.empty()) {
492  // get the pointer to the energy tensor, dimension is batch x 1
493  float *energy = outputs[0].flat<float>().data();
494 
495  for (const int &i : tracksterIndices) {
496  tracksters[i].setRegressedEnergy(*(energy++));
497  }
498  }
499 
500  // store id probabilities per trackster (8)
501  if (!eidOutputNameId_.empty()) {
502  // get the pointer to the id probability tensor, dimension is batch x id_probabilities.size()
503  int probsIdx = eidOutputNameEnergy_.empty() ? 0 : 1;
504  float *probs = outputs[probsIdx].flat<float>().data();
505 
506  for (const int &i : tracksterIndices) {
507  tracksters[i].setProbabilities(probs);
508  probs += tracksters[i].id_probabilities().size();
509  }
510  }
511 }
std::vector< NamedTensor > NamedTensorList
Definition: TensorFlow.h:31
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:184
static std::string const input
Definition: EdmProvDump.cc:50
std::vector< float > features(const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:281
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
string inputList
Definition: crabTemplate.py:6
double energy() const
cluster energy
Definition: CaloCluster.h:149
std::vector< unsigned int > & vertices()
Definition: Trackster.h:57
std::vector< float > & vertex_multiplicity()
Definition: Trackster.h:58
double b
Definition: hdecay.h:120
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
double a
Definition: hdecay.h:121
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381

◆ fillPSetDescription()

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

Definition at line 859 of file PatternRecognitionbyCLUE3D.cc.

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

859  {
860  iDesc.add<int>("algo_verbosity", 0);
861  iDesc.add<std::vector<double>>("criticalDensity", {4, 4, 4})->setComment("in GeV");
862  iDesc.add<std::vector<double>>("criticalSelfDensity", {0.15, 0.15, 0.15} /* roughly 1/(densitySiblingLayers+1) */)
863  ->setComment("Minimum ratio of self_energy/local_density to become a seed.");
864  iDesc.add<std::vector<int>>("densitySiblingLayers", {3, 3, 3})
865  ->setComment(
866  "inclusive, layers to consider while computing local density and searching for nearestHigher higher");
867  iDesc.add<std::vector<double>>("densityEtaPhiDistanceSqr", {0.0008, 0.0008, 0.0008})
868  ->setComment("in eta,phi space, distance to consider for local density");
869  iDesc.add<std::vector<double>>("densityXYDistanceSqr", {3.24, 3.24, 3.24})
870  ->setComment("in cm, distance on the transverse plane to consider for local density");
871  iDesc.add<std::vector<double>>("kernelDensityFactor", {0.2, 0.2, 0.2})
872  ->setComment("Kernel factor to be applied to other LC while computing the local density");
873  iDesc.add<bool>("densityOnSameLayer", false);
874  iDesc.add<bool>("nearestHigherOnSameLayer", false)
875  ->setComment("Allow the nearestHigher to be located on the same layer");
876  iDesc.add<bool>("useAbsoluteProjectiveScale", true)
877  ->setComment("Express all cuts in terms of r/z*z_0{,phi} projective variables");
878  iDesc.add<bool>("useClusterDimensionXY", false)
879  ->setComment(
880  "Boolean. If true use the estimated cluster radius to determine the cluster compatibility while computing "
881  "the local density");
882  iDesc.add<bool>("rescaleDensityByZ", false)
883  ->setComment(
884  "Rescale local density by the extension of the Z 'volume' explored. The transvere dimension is, at present, "
885  "fixed and factored out.");
886  iDesc.add<std::vector<double>>("criticalEtaPhiDistance", {0.025, 0.025, 0.025})
887  ->setComment("Minimal distance in eta,phi space from nearestHigher to become a seed");
888  iDesc.add<std::vector<double>>("criticalXYDistance", {1.8, 1.8, 1.8})
889  ->setComment("Minimal distance in cm on the XY plane from nearestHigher to become a seed");
890  iDesc.add<std::vector<int>>("criticalZDistanceLyr", {5, 5, 5})
891  ->setComment("Minimal distance in layers along the Z axis from nearestHigher to become a seed");
892  iDesc.add<std::vector<double>>("outlierMultiplier", {2, 2, 2})
893  ->setComment("Minimal distance in transverse space from nearestHigher to become an outlier");
894  iDesc.add<std::vector<int>>("minNumLayerCluster", {2, 2, 2})->setComment("Not Inclusive");
895  iDesc.add<bool>("doPidCut", false);
896  iDesc.add<double>("cutHadProb", 0.5);
897  iDesc.add<std::string>("eid_input_name", "input");
898  iDesc.add<std::string>("eid_output_name_energy", "output/regressed_energy");
899  iDesc.add<std::string>("eid_output_name_id", "output/id_probabilities");
900  iDesc.add<double>("eid_min_cluster_energy", 1.);
901  iDesc.add<int>("eid_n_layers", 50);
902  iDesc.add<int>("eid_n_clusters", 10);
903  iDesc.add<bool>("computeLocalTime", false);
904 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ findAndAssignTracksters()

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

Definition at line 786 of file PatternRecognitionbyCLUE3D.cc.

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

787  {
788  unsigned int nTracksters = 0;
789  std::vector<std::pair<int, int>> localStack;
790  const auto &critical_transverse_distance =
792  // find cluster seeds and outlier
793  for (unsigned int layer = 0; layer < 2 * rhtools_.lastLayer(); layer++) {
794  auto &clustersOnLayer = clusters_[layer];
795  unsigned int numberOfClusters = clustersOnLayer.x.size();
796  for (unsigned int i = 0; i < numberOfClusters; i++) {
797  // initialize clusterIndex
798  clustersOnLayer.clusterIndex[i] = -1;
799  auto algoId = clustersOnLayer.algoId[i];
800  bool isSeed = (clustersOnLayer.delta[i].first > critical_transverse_distance[algoId] ||
801  clustersOnLayer.delta[i].second > criticalZDistanceLyr_[algoId]) &&
802  (clustersOnLayer.rho[i] >= criticalDensity_[algoId]) &&
803  (clustersOnLayer.energy[i] / clustersOnLayer.rho[i] > criticalSelfDensity_[algoId]);
804  if (!clustersOnLayer.isSilicon[i]) {
805  isSeed = (clustersOnLayer.delta[i].first > clustersOnLayer.radius[i] ||
806  clustersOnLayer.delta[i].second > criticalZDistanceLyr_[algoId]) &&
807  (clustersOnLayer.rho[i] >= criticalDensity_[algoId]) &&
808  (clustersOnLayer.energy[i] / clustersOnLayer.rho[i] > criticalSelfDensity_[algoId]);
809  }
810  bool isOutlier =
811  (clustersOnLayer.delta[i].first > outlierMultiplier_[algoId] * critical_transverse_distance[algoId]) &&
812  (clustersOnLayer.rho[i] < criticalDensity_[algoId]);
813  if (isSeed) {
815  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
816  << "Found seed on Layer " << layer << " SOAidx: " << i << " assigned ClusterIdx: " << nTracksters;
817  }
818  clustersOnLayer.clusterIndex[i] = nTracksters++;
819  tracksterSeedAlgoId_.push_back(algoId);
820  clustersOnLayer.isSeed[i] = true;
821  localStack.emplace_back(layer, i);
822  } else if (!isOutlier) {
823  auto [lyrIdx, soaIdx] = clustersOnLayer.nearestHigher[i];
825  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
826  << "Found follower on Layer " << layer << " SOAidx: " << i << " attached to cluster on layer: " << lyrIdx
827  << " SOAidx: " << soaIdx;
828  }
829  if (lyrIdx >= 0)
830  clusters_[lyrIdx].followers[soaIdx].emplace_back(layer, i);
831  } else {
833  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
834  << "Found Outlier on Layer " << layer << " SOAidx: " << i << " with rho: " << clustersOnLayer.rho[i]
835  << " and delta: " << clustersOnLayer.delta[i].first << ", " << clustersOnLayer.delta[i].second;
836  }
837  }
838  }
839  }
840 
841  // Propagate cluster index
842  while (!localStack.empty()) {
843  auto [lyrIdx, soaIdx] = localStack.back();
844  auto &thisSeed = clusters_[lyrIdx].followers[soaIdx];
845  localStack.pop_back();
846 
847  // loop over followers
848  for (auto [follower_lyrIdx, follower_soaIdx] : thisSeed) {
849  // pass id to a follower
850  clusters_[follower_lyrIdx].clusterIndex[follower_soaIdx] = clusters_[lyrIdx].clusterIndex[soaIdx];
851  // push this follower to localStack
852  localStack.emplace_back(follower_lyrIdx, follower_soaIdx);
853  }
854  }
855  return nTracksters;
856 }
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 154 of file PatternRecognitionbyCLUE3D.cc.

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

157  {
158  // Protect from events with no seeding regions
159  if (input.regions.empty())
160  return;
161 
162  const int eventNumber = input.ev.eventAuxiliary().event();
164  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "New Event";
165  }
166 
167  edm::EventSetup const &es = input.es;
170 
171  // Assume identical Z-positioning between positive and negative sides.
172  // Also, layers inside the HGCAL geometry start from 1.
173  for (unsigned int i = 0; i < rhtools_.lastLayer(); ++i) {
174  layersPosZ_.push_back(rhtools_.getPositionLayer(i + 1).z());
176  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Layer " << i << " located at Z: " << layersPosZ_.back();
177  }
178  }
179 
180  clusters_.clear();
181  tracksterSeedAlgoId_.clear();
182 
183  clusters_.resize(2 * rhtools_.lastLayer(false));
184  std::vector<std::pair<int, int>> layerIdx2layerandSoa; //used everywhere also to propagate cluster masking
185 
186  layerIdx2layerandSoa.reserve(input.layerClusters.size());
187  unsigned int layerIdx = 0;
188  for (auto const &lc : input.layerClusters) {
189  if (input.mask[layerIdx] == 0.) {
191  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Skipping masked cluster: " << layerIdx;
192  }
193  layerIdx2layerandSoa.emplace_back(-1, -1);
194  layerIdx++;
195  continue;
196  }
197  const auto firstHitDetId = lc.hitsAndFractions()[0].first;
198  int layer = rhtools_.getLayerWithOffset(firstHitDetId) - 1 +
199  rhtools_.lastLayer(false) * ((rhtools_.zside(firstHitDetId) + 1) >> 1);
200  assert(layer >= 0);
201  auto detId = lc.hitsAndFractions()[0].first;
202  int layerClusterIndexInLayer = clusters_[layer].x.size();
203  layerIdx2layerandSoa.emplace_back(layer, layerClusterIndexInLayer);
204  float sum_x = 0.;
205  float sum_y = 0.;
206  float sum_sqr_x = 0.;
207  float sum_sqr_y = 0.;
208  float ref_x = lc.x();
209  float ref_y = lc.y();
210  float invClsize = 1. / lc.hitsAndFractions().size();
211  for (auto const &hitsAndFractions : lc.hitsAndFractions()) {
212  auto const &point = rhtools_.getPosition(hitsAndFractions.first);
213  sum_x += point.x() - ref_x;
214  sum_sqr_x += (point.x() - ref_x) * (point.x() - ref_x);
215  sum_y += point.y() - ref_y;
216  sum_sqr_y += (point.y() - ref_y) * (point.y() - ref_y);
217  }
218  // The variance of X for X uniform in circle of radius R is R^2/4,
219  // therefore we multiply the sqrt(var) by 2 to have a rough estimate of the
220  // radius. On the other hand, while averaging the x and y radius, we would
221  // end up dividing by 2. Hence we omit the value here and in the average
222  // below, too.
223  float radius_x = sqrt((sum_sqr_x - (sum_x * sum_x) * invClsize) * invClsize);
224  float radius_y = sqrt((sum_sqr_y - (sum_y * sum_y) * invClsize) * invClsize);
226  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
227  << "cluster rx: " << std::setw(5) << radius_x << ", ry: " << std::setw(5) << radius_y
228  << ", r: " << std::setw(5) << (radius_x + radius_y) << ", cells: " << std::setw(4)
229  << lc.hitsAndFractions().size();
230  }
231 
232  // The case of single cell layer clusters has to be handled differently.
233 
234  if (invClsize == 1.) {
235  // Silicon case
236  if (rhtools_.isSilicon(detId)) {
237  radius_x = radius_y = rhtools_.getRadiusToSide(detId);
239  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Single cell cluster in silicon, rx: " << std::setw(5)
240  << radius_x << ", ry: " << std::setw(5) << radius_y;
241  }
242  } else {
243  auto const &point = rhtools_.getPosition(detId);
244  auto const &eta_phi_window = rhtools_.getScintDEtaDPhi(detId);
245  radius_x = radius_y = point.perp() * eta_phi_window.second;
247  edm::LogVerbatim("PatternRecognitionbyCLUE3D")
248  << "Single cell cluster in scintillator. rx: " << std::setw(5) << radius_x << ", ry: " << std::setw(5)
249  << radius_y << ", eta-span: " << std::setw(5) << eta_phi_window.first << ", phi-span: " << std::setw(5)
250  << eta_phi_window.second;
251  }
252  }
253  }
254 
255  // Maybe check if these vectors can be reserved beforehands
256  clusters_[layer].x.emplace_back(lc.x());
257  clusters_[layer].y.emplace_back(lc.y());
258  clusters_[layer].z.emplace_back(lc.z());
259  clusters_[layer].r_over_absz.emplace_back(sqrt(lc.x() * lc.x() + lc.y() * lc.y()) / std::abs(lc.z()));
260  clusters_[layer].radius.emplace_back(radius_x + radius_y);
261  clusters_[layer].eta.emplace_back(lc.eta());
262  clusters_[layer].phi.emplace_back(lc.phi());
263  clusters_[layer].cells.push_back(lc.hitsAndFractions().size());
264  clusters_[layer].algoId.push_back(lc.algo() - reco::CaloCluster::hgcal_em);
265  clusters_[layer].isSilicon.push_back(rhtools_.isSilicon(detId));
266  clusters_[layer].energy.emplace_back(lc.energy());
267  clusters_[layer].isSeed.push_back(false);
268  clusters_[layer].clusterIndex.emplace_back(-1);
269  clusters_[layer].layerClusterOriginalIdx.emplace_back(layerIdx++);
270  clusters_[layer].nearestHigher.emplace_back(-1, -1);
271  clusters_[layer].rho.emplace_back(0.f);
272  clusters_[layer].z_extension.emplace_back(0.f);
273  clusters_[layer].delta.emplace_back(
275  }
276  for (unsigned int layer = 0; layer < clusters_.size(); layer++) {
277  clusters_[layer].followers.resize(clusters_[layer].x.size());
278  }
279 
280  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
281  int maxLayer = 2 * lastLayerPerSide - 1;
282  std::vector<int> numberOfClustersPerLayer(maxLayer, 0);
283  for (int i = 0; i <= maxLayer; i++) {
284  calculateLocalDensity(input.tiles, i, layerIdx2layerandSoa);
285  }
286  for (int i = 0; i <= maxLayer; i++) {
287  calculateDistanceToHigher(input.tiles, i, layerIdx2layerandSoa);
288  }
289 
290  auto nTracksters = findAndAssignTracksters(input.tiles, layerIdx2layerandSoa);
292  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Reconstructed " << nTracksters << " tracksters" << std::endl;
293  dumpClusters(input.tiles, layerIdx2layerandSoa, eventNumber);
294  }
295 
296  // Build Trackster
297  result.resize(nTracksters);
298 
299  for (unsigned int layer = 0; layer < clusters_.size(); ++layer) {
300  const auto &thisLayer = clusters_[layer];
302  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Examining Layer: " << layer;
303  }
304  for (unsigned int lc = 0; lc < thisLayer.x.size(); ++lc) {
306  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Trackster " << thisLayer.clusterIndex[lc];
307  }
308  if (thisLayer.clusterIndex[lc] >= 0) {
310  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << " adding lcIdx: " << thisLayer.layerClusterOriginalIdx[lc];
311  }
312  result[thisLayer.clusterIndex[lc]].vertices().push_back(thisLayer.layerClusterOriginalIdx[lc]);
313  result[thisLayer.clusterIndex[lc]].vertex_multiplicity().push_back(1);
314  // loop over followers
315  for (auto [follower_lyrIdx, follower_soaIdx] : thisLayer.followers[lc]) {
316  std::array<unsigned int, 2> edge = {
317  {(unsigned int)thisLayer.layerClusterOriginalIdx[lc],
318  (unsigned int)clusters_[follower_lyrIdx].layerClusterOriginalIdx[follower_soaIdx]}};
319  result[thisLayer.clusterIndex[lc]].edges().push_back(edge);
320  }
321  }
322  }
323  }
324  size_t tracksterIndex = 0;
325  result.erase(std::remove_if(std::begin(result),
326  std::end(result),
327  [&](auto const &v) {
328  return static_cast<int>(v.vertices().size()) <
329  minNumLayerCluster_[tracksterSeedAlgoId_[tracksterIndex++]];
330  }),
331  result.end());
332  if (doPidCut_) {
333  energyRegressionAndID(input.layerClusters, input.tfSession, result);
334  result.erase(std::remove_if(std::begin(result),
335  std::end(result),
336  [&](auto const &v) {
337  auto const &hadProb =
340  return hadProb >= cutHadProb_;
341  }),
342  result.end());
343  }
344  result.shrink_to_fit();
345 
347  input.layerClusters,
348  input.layerClustersTime,
351 
353  for (auto const &t : result) {
354  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Barycenter: " << t.barycenter();
355  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "LCs: " << t.vertices().size();
356  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Energy: " << t.raw_energy();
357  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << "Regressed: " << t.regressed_energy();
358  }
359  }
360 
361  // Dump Tracksters information
363  dumpTracksters(layerIdx2layerandSoa, eventNumber, result);
364  }
365 
366  // Reset internal clusters_ structure of array for next event
367  reset();
369  edm::LogVerbatim("PatternRecognitionbyCLUE3D") << std::endl;
370  }
371 }
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
void energyRegressionAndID(const std::vector< reco::CaloCluster > &layerClusters, const tensorflow::Session *, std::vector< Trackster > &result)
T z() const
Definition: PV3DBase.h:61
std::float_t getRadiusToSide(const DetId &) const
Definition: RecHitTools.cc:246
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 >> &)
void assignPCAtoTracksters(std::vector< Trackster > &, const std::vector< reco::CaloCluster > &, const edm::ValueMap< std::pair< float, float >> &, double, bool computeLocalTime=false, bool energyWeight=true)
*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 96 of file PatternRecognitionbyCLUE3D.h.

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

96  {
97  for (auto& c : clusters_) {
98  c.clear();
99  c.shrink_to_fit();
100  }
101  }
std::vector< ClustersOnLayer > clusters_

Member Data Documentation

◆ caloGeomToken_

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

Definition at line 117 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 143 of file PatternRecognitionbyCLUE3D.h.

◆ criticalDensity_

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

Definition at line 118 of file PatternRecognitionbyCLUE3D.h.

◆ criticalEtaPhiDistance_

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

Definition at line 129 of file PatternRecognitionbyCLUE3D.h.

◆ criticalSelfDensity_

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

Definition at line 119 of file PatternRecognitionbyCLUE3D.h.

◆ criticalXYDistance_

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

Definition at line 130 of file PatternRecognitionbyCLUE3D.h.

◆ criticalZDistanceLyr_

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

Definition at line 131 of file PatternRecognitionbyCLUE3D.h.

◆ cutHadProb_

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

Definition at line 135 of file PatternRecognitionbyCLUE3D.h.

◆ densityEtaPhiDistanceSqr_

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

Definition at line 121 of file PatternRecognitionbyCLUE3D.h.

◆ densityOnSameLayer_

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

Definition at line 124 of file PatternRecognitionbyCLUE3D.h.

◆ densitySiblingLayers_

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

Definition at line 120 of file PatternRecognitionbyCLUE3D.h.

◆ densityXYDistanceSqr_

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

Definition at line 122 of file PatternRecognitionbyCLUE3D.h.

◆ doPidCut_

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

Definition at line 134 of file PatternRecognitionbyCLUE3D.h.

◆ eidInputName_

template<typename TILES >
const std::string ticl::PatternRecognitionbyCLUE3D< TILES >::eidInputName_
private

Definition at line 137 of file PatternRecognitionbyCLUE3D.h.

◆ eidMinClusterEnergy_

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

Definition at line 140 of file PatternRecognitionbyCLUE3D.h.

◆ eidNClusters_

template<typename TILES >
const int ticl::PatternRecognitionbyCLUE3D< TILES >::eidNClusters_
private

Definition at line 142 of file PatternRecognitionbyCLUE3D.h.

◆ eidNFeatures_

template<typename TILES >
const int ticl::PatternRecognitionbyCLUE3D< TILES >::eidNFeatures_ = 3
staticprivate

Definition at line 148 of file PatternRecognitionbyCLUE3D.h.

◆ eidNLayers_

template<typename TILES >
const int ticl::PatternRecognitionbyCLUE3D< TILES >::eidNLayers_
private

Definition at line 141 of file PatternRecognitionbyCLUE3D.h.

◆ eidOutputNameEnergy_

template<typename TILES >
const std::string ticl::PatternRecognitionbyCLUE3D< TILES >::eidOutputNameEnergy_
private

Definition at line 138 of file PatternRecognitionbyCLUE3D.h.

◆ eidOutputNameId_

template<typename TILES >
const std::string ticl::PatternRecognitionbyCLUE3D< TILES >::eidOutputNameId_
private

Definition at line 139 of file PatternRecognitionbyCLUE3D.h.

◆ eidSession_

template<typename TILES >
tensorflow::Session* ticl::PatternRecognitionbyCLUE3D< TILES >::eidSession_
private

Definition at line 146 of file PatternRecognitionbyCLUE3D.h.

◆ filter_on_categories_

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

Definition at line 136 of file PatternRecognitionbyCLUE3D.h.

◆ kernelDensityFactor_

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

Definition at line 123 of file PatternRecognitionbyCLUE3D.h.

◆ layersPosZ_

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

Definition at line 114 of file PatternRecognitionbyCLUE3D.h.

◆ minNumLayerCluster_

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

Definition at line 133 of file PatternRecognitionbyCLUE3D.h.

◆ nearestHigherOnSameLayer_

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

Definition at line 125 of file PatternRecognitionbyCLUE3D.h.

◆ outlierMultiplier_

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

Definition at line 132 of file PatternRecognitionbyCLUE3D.h.

◆ rescaleDensityByZ_

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

Definition at line 128 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 115 of file PatternRecognitionbyCLUE3D.h.

◆ useAbsoluteProjectiveScale_

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

Definition at line 126 of file PatternRecognitionbyCLUE3D.h.

◆ useClusterDimensionXY_

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

Definition at line 127 of file PatternRecognitionbyCLUE3D.h.