CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 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,
CaloGeometryRecord
caloGeomToken_
 
std::vector< ClustersOnLayerclusters_
 
const double criticalDensity_
 
const double criticalEtaPhiDistance_
 
const double densityEtaPhiDistanceSqr_
 
const double densityOnSameLayer_
 
const int densitySiblingLayers_
 
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 int minNumLayerCluster_
 
const double outlierMultiplier_
 
hgcal::RecHitTools rhtools_
 

Static Private Attributes

static const int eidNFeatures_ = 3
 

Additional Inherited Members

- Public Types inherited from ticl::PatternRecognitionAlgoBaseT< TILES >
enum  VerbosityLevel {
  None = 0, Basic, Advanced, Expert,
  Guru
}
 
- 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

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

Definition at line 21 of file PatternRecognitionbyCLUE3D.cc.

24  criticalDensity_(conf.getParameter<double>("criticalDensity")),
25  densitySiblingLayers_(conf.getParameter<int>("densitySiblingLayers")),
26  densityEtaPhiDistanceSqr_(conf.getParameter<double>("densityEtaPhiDistanceSqr")),
27  densityOnSameLayer_(conf.getParameter<bool>("densityOnSameLayer")),
28  criticalEtaPhiDistance_(conf.getParameter<double>("criticalEtaPhiDistance")),
29  outlierMultiplier_(conf.getParameter<double>("outlierMultiplier")),
30  minNumLayerCluster_(conf.getParameter<int>("minNumLayerCluster")),
31  eidInputName_(conf.getParameter<std::string>("eid_input_name")),
32  eidOutputNameEnergy_(conf.getParameter<std::string>("eid_output_name_energy")),
33  eidOutputNameId_(conf.getParameter<std::string>("eid_output_name_id")),
34  eidMinClusterEnergy_(conf.getParameter<double>("eid_min_cluster_energy")),
35  eidNLayers_(conf.getParameter<int>("eid_n_layers")),
36  eidNClusters_(conf.getParameter<int>("eid_n_clusters")){};
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
template<typename TILES >
ticl::PatternRecognitionbyCLUE3D< TILES >::~PatternRecognitionbyCLUE3D ( )
overridedefault

Member Function Documentation

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

Definition at line 548 of file PatternRecognitionbyCLUE3D.cc.

References reco::deltaR2(), mps_fire::i, SiStripPI::max, SiStripPI::min, HLT_FULL_cff::minLayer, L1TMuonDQMOffline_cfi::nEtaBins, ecaldqm::binning::nPhiBins, hltrates_dqm_sourceclient-live_cfg::offset, and mathSSE::sqrt().

549  {
550  constexpr int nEtaBin = TILES::constants_type_t::nEtaBins;
551  constexpr int nPhiBin = TILES::constants_type_t::nPhiBins;
552  auto &clustersOnLayer = clusters_[layerId];
553  unsigned int numberOfClusters = clustersOnLayer.x.size();
554 
555  for (unsigned int i = 0; i < numberOfClusters; i++) {
557  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
558  << "Starting searching nearestHigher on " << layerId << " with rho: " << clustersOnLayer.rho[i]
559  << " at eta, phi: " << tiles[layerId].etaBin(clustersOnLayer.eta[i]) << ", "
560  << tiles[layerId].etaBin(clustersOnLayer.phi[i]);
561  }
562  // We need to partition the two sides of the HGCAL detector
563  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
564  int minLayer = 0;
565  int maxLayer = 2 * lastLayerPerSide - 1;
566  if (layerId < lastLayerPerSide) {
567  minLayer = std::max(layerId - densitySiblingLayers_, minLayer);
568  maxLayer = std::min(layerId + densitySiblingLayers_, lastLayerPerSide - 1);
569  } else {
570  minLayer = std::max(layerId - densitySiblingLayers_, lastLayerPerSide + 1);
571  maxLayer = std::min(layerId + densitySiblingLayers_, maxLayer);
572  }
573  float maxDelta = std::numeric_limits<float>::max();
574  float i_delta = maxDelta;
575  std::pair<int, int> i_nearestHigher(-1, -1);
576  for (int currentLayer = minLayer; currentLayer <= maxLayer; currentLayer++) {
577  const auto &tileOnLayer = tiles[currentLayer];
578  int etaWindow = 3;
579  int phiWindow = 3;
580  int etaBinMin = std::max(tileOnLayer.etaBin(clustersOnLayer.eta[i]) - etaWindow, 0);
581  int etaBinMax = std::min(tileOnLayer.etaBin(clustersOnLayer.eta[i]) + etaWindow, nEtaBin);
582  int phiBinMin = tileOnLayer.phiBin(clustersOnLayer.phi[i]) - phiWindow;
583  int phiBinMax = tileOnLayer.phiBin(clustersOnLayer.phi[i]) + phiWindow;
584  for (int ieta = etaBinMin; ieta <= etaBinMax; ++ieta) {
585  auto offset = ieta * nPhiBin;
586  for (int iphi_it = phiBinMin; iphi_it <= phiBinMax; ++iphi_it) {
587  int iphi = ((iphi_it % nPhiBin + nPhiBin) % nPhiBin);
589  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
590  << "Searching nearestHigher on " << currentLayer << " eta, phi: " << ieta << ", " << iphi_it;
591  }
592  for (auto otherClusterIdx : tileOnLayer[offset + iphi]) {
593  auto const &layerandSoa = layerIdx2layerandSoa[otherClusterIdx];
594  // Skip masked layer clusters
595  if ((layerandSoa.first == -1) && (layerandSoa.second == -1))
596  continue;
597  auto const &clustersOnOtherLayer = clusters_[layerandSoa.first];
598  float dist = reco::deltaR2(clustersOnLayer.eta[i],
599  clustersOnLayer.phi[i],
600  clustersOnOtherLayer.eta[layerandSoa.second],
601  clustersOnOtherLayer.phi[layerandSoa.second]);
602  bool foundHigher = (clustersOnOtherLayer.rho[layerandSoa.second] > clustersOnLayer.rho[i]) ||
603  (clustersOnOtherLayer.rho[layerandSoa.second] == clustersOnLayer.rho[i] &&
604  clustersOnOtherLayer.layerClusterOriginalIdx[layerandSoa.second] >
605  clustersOnLayer.layerClusterOriginalIdx[i]);
607  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
608  << "Searching nearestHigher on " << currentLayer
609  << " with rho: " << clustersOnOtherLayer.rho[layerandSoa.second]
610  << " on layerIdxInSOA: " << layerandSoa.first << ", " << layerandSoa.second
611  << " with distance: " << sqrt(dist) << " foundHigher: " << foundHigher;
612  }
613  if (foundHigher && dist <= i_delta) {
614  // update i_delta
615  i_delta = dist;
616  // update i_nearestHigher
617  i_nearestHigher = layerandSoa;
618  }
619  } // End of loop on clusters
620  } // End of loop on phi bins
621  } // End of loop on eta bins
622  } // End of loop on layers
623 
624  bool foundNearestHigherInEtaPhiCylinder = (i_delta != maxDelta);
626  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
627  << "i_delta: " << sqrt(i_delta) << " passed: " << foundNearestHigherInEtaPhiCylinder << " "
628  << i_nearestHigher.first << " " << i_nearestHigher.second;
629  }
630  if (foundNearestHigherInEtaPhiCylinder) {
631  clustersOnLayer.delta[i] = sqrt(i_delta);
632  clustersOnLayer.nearestHigher[i] = i_nearestHigher;
633  } else {
634  // otherwise delta is guaranteed to be larger outlierDeltaFactor_*delta_c
635  // we can safely maximize delta to be maxDelta
636  clustersOnLayer.delta[i] = maxDelta;
637  clustersOnLayer.nearestHigher[i] = {-1, -1};
638  }
639  } // End of loop on clusters
640 }
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:79
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
template<typename TILES >
void PatternRecognitionbyCLUE3D::calculateLocalDensity ( const TILES &  tiles,
const int  layerId,
const std::vector< std::pair< int, int >> &  layerIdx2layerandSoa 
)
private

Definition at line 429 of file PatternRecognitionbyCLUE3D.cc.

References CommonMethods::delta(), reco::deltaR2(), mps_fire::i, SiStripPI::max, SiStripPI::min, HLT_FULL_cff::minLayer, L1TMuonDQMOffline_cfi::nEtaBins, ecaldqm::binning::nPhiBins, and hltrates_dqm_sourceclient-live_cfg::offset.

430  {
431  constexpr int nEtaBin = TILES::constants_type_t::nEtaBins;
432  constexpr int nPhiBin = TILES::constants_type_t::nPhiBins;
433  auto &clustersOnLayer = clusters_[layerId];
434  unsigned int numberOfClusters = clustersOnLayer.x.size();
435 
436  auto isReachable = [&](float x1, float x2, float y1, float y2, float delta_sqr) -> bool {
438  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
439  << ((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2)) << " vs " << delta_sqr << "["
440  << ((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2) < delta_sqr) << "]\n";
441  }
442  return (x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2) < delta_sqr;
443  };
444 
445  for (unsigned int i = 0; i < numberOfClusters; i++) {
446  // We need to partition the two sides of the HGCAL detector
447  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
448  int minLayer = 0;
449  int maxLayer = 2 * lastLayerPerSide - 1;
450  if (layerId < lastLayerPerSide) {
451  minLayer = std::max(layerId - densitySiblingLayers_, minLayer);
452  maxLayer = std::min(layerId + densitySiblingLayers_, lastLayerPerSide - 1);
453  } else {
454  minLayer = std::max(layerId - densitySiblingLayers_, lastLayerPerSide);
455  maxLayer = std::min(layerId + densitySiblingLayers_, maxLayer);
456  }
457  for (int currentLayer = minLayer; currentLayer <= maxLayer; currentLayer++) {
459  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "RefLayer: " << layerId << " SoaIDX: " << i;
460  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "NextLayer: " << currentLayer;
461  }
462  const auto &tileOnLayer = tiles[currentLayer];
463  bool onSameLayer = (currentLayer == layerId);
465  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "onSameLayer: " << onSameLayer;
466  }
467  const int etaWindow = 2;
468  const int phiWindow = 2;
469  int etaBinMin = std::max(tileOnLayer.etaBin(clustersOnLayer.eta[i]) - etaWindow, 0);
470  int etaBinMax = std::min(tileOnLayer.etaBin(clustersOnLayer.eta[i]) + etaWindow, nEtaBin);
471  int phiBinMin = tileOnLayer.phiBin(clustersOnLayer.phi[i]) - phiWindow;
472  int phiBinMax = tileOnLayer.phiBin(clustersOnLayer.phi[i]) + phiWindow;
474  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "eta: " << clustersOnLayer.eta[i];
475  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "phi: " << clustersOnLayer.phi[i];
476  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "etaBinMin: " << etaBinMin << ", etaBinMax: " << etaBinMax;
477  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "phiBinMin: " << phiBinMin << ", phiBinMax: " << phiBinMax;
478  }
479  for (int ieta = etaBinMin; ieta <= etaBinMax; ++ieta) {
480  auto offset = ieta * nPhiBin;
482  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "offset: " << offset;
483  }
484  for (int iphi_it = phiBinMin; iphi_it <= phiBinMax; ++iphi_it) {
485  int iphi = ((iphi_it % nPhiBin + nPhiBin) % nPhiBin);
487  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "iphi: " << iphi;
488  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
489  << "Entries in tileBin: " << tileOnLayer[offset + iphi].size();
490  }
491  for (auto otherClusterIdx : tileOnLayer[offset + iphi]) {
492  auto const &layerandSoa = layerIdx2layerandSoa[otherClusterIdx];
493  // Skip masked layer clusters
494  if ((layerandSoa.first == -1) && (layerandSoa.second == -1)) {
496  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Skipping masked layerIdx " << otherClusterIdx;
497  }
498  continue;
499  }
500  auto const &clustersLayer = clusters_[layerandSoa.first];
502  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
503  << "OtherLayer: " << layerandSoa.first << " SoaIDX: " << layerandSoa.second;
504  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "OtherEta: " << clustersLayer.eta[layerandSoa.second];
505  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "OtherPhi: " << clustersLayer.phi[layerandSoa.second];
506  }
507  // extend by 26 mm, roughly 2 cells, more wrt sum of radii
508  float delta = clustersOnLayer.radius[i] + clustersLayer.radius[layerandSoa.second] + 2.6f;
509  if (densityOnSameLayer_ && onSameLayer) {
510  if (isReachable(clustersOnLayer.x[i],
511  clustersLayer.x[layerandSoa.second],
512  clustersOnLayer.y[i],
513  clustersLayer.y[layerandSoa.second],
514  delta * delta)) {
515  clustersOnLayer.rho[i] += (clustersOnLayer.layerClusterOriginalIdx[i] == otherClusterIdx ? 1.f : 0.2f) *
516  clustersLayer.energy[layerandSoa.second];
517  }
518  } else {
520  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Distance: "
521  << reco::deltaR2(clustersOnLayer.eta[i],
522  clustersOnLayer.phi[i],
523  clustersLayer.eta[layerandSoa.second],
524  clustersLayer.phi[layerandSoa.second]);
525  }
526  if (reco::deltaR2(clustersOnLayer.eta[i],
527  clustersOnLayer.phi[i],
528  clustersLayer.eta[layerandSoa.second],
529  clustersLayer.phi[layerandSoa.second]) < densityEtaPhiDistanceSqr_) {
530  auto energyToAdd = (clustersOnLayer.layerClusterOriginalIdx[i] == otherClusterIdx ? 1.f : 0.5f) *
531  clustersLayer.energy[layerandSoa.second];
532  clustersOnLayer.rho[i] += energyToAdd;
533  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
534  << "Adding " << energyToAdd << " partial " << clustersOnLayer.rho[i];
535  }
536  }
537  } // end of loop on possible compatible clusters
538  } // end of loop over phi-bin region
539  } // end of loop over eta-bin region
540  } // end of loop on the sibling layers
541  } // end of loop over clusters on this layer
543  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << std::endl;
544  }
545 }
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:79
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
template<typename TILES >
void PatternRecognitionbyCLUE3D::dumpClusters ( const std::vector< std::pair< int, int >> &  layerIdx2layerandSoa,
const int  eventNumber 
) const
private

Definition at line 93 of file PatternRecognitionbyCLUE3D.cc.

References phase1PixelTopology::layer, pileupDistInMC::num, and findQualityFiles::v.

94  {
96  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "[evt, layer, x, y, eta, phi, cells, energy, radius, rho, delta, "
97  "isSeed, clusterIdx, layerClusterOriginalIdx";
98  }
99 
100  for (unsigned int layer = 0; layer < clusters_.size(); layer++) {
102  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "On Layer " << layer;
103  }
104  auto const &thisLayer = clusters_[layer];
105  int num = 0;
106  for (auto v : thisLayer.x) {
108  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
109  << "ClusterInfo: " << eventNumber << ", " << layer << ", " << v << ", " << thisLayer.y[num] << ", "
110  << thisLayer.eta[num] << ", " << thisLayer.phi[num] << ", " << thisLayer.cells[num] << ", "
111  << thisLayer.energy[num] << ", " << thisLayer.radius[num] << ", " << thisLayer.rho[num] << ", "
112  << thisLayer.delta[num] << ", " << thisLayer.isSeed[num] << ", " << thisLayer.clusterIndex[num] << ", "
113  << thisLayer.layerClusterOriginalIdx[num];
114  }
115  ++num;
116  }
117  }
118  for (unsigned int lcIdx = 0; lcIdx < layerIdx2layerandSoa.size(); lcIdx++) {
119  auto const &layerandSoa = layerIdx2layerandSoa[lcIdx];
120  // Skip masked layer clusters
121  if ((layerandSoa.first == -1) && (layerandSoa.second == -1))
122  continue;
124  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
125  << "lcIdx: " << lcIdx << " on Layer: " << layerandSoa.first << " SOA: " << layerandSoa.second;
126  }
127  }
128 }
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_
constexpr std::array< uint8_t, layerIndexSize > layer
template<typename TILES >
void PatternRecognitionbyCLUE3D::dumpTiles ( const TILES &  tiles) const
private

Definition at line 39 of file PatternRecognitionbyCLUE3D.cc.

References relativeConstraints::empty, phase1PixelTopology::layer, L1TMuonDQMOffline_cfi::nEtaBins, ecaldqm::binning::nPhiBins, and hltrates_dqm_sourceclient-live_cfg::offset.

39  {
40  constexpr int nEtaBin = TILES::constants_type_t::nEtaBins;
41  constexpr int nPhiBin = TILES::constants_type_t::nPhiBins;
42  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
43  int maxLayer = 2 * lastLayerPerSide - 1;
44  for (int layer = 0; layer <= maxLayer; layer++) {
45  for (int ieta = 0; ieta < nEtaBin; ieta++) {
46  auto offset = ieta * nPhiBin;
47  for (int phi = 0; phi < nPhiBin; phi++) {
48  int iphi = ((phi % nPhiBin + nPhiBin) % nPhiBin);
49  if (!tiles[layer][offset + iphi].empty()) {
50  if (this->algo_verbosity_ > this->Advanced) {
51  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Layer: " << layer << " ieta: " << ieta << " phi: " << phi
52  << " " << tiles[layer][offset + iphi].size();
53  }
54  }
55  }
56  }
57  }
58 }
Log< level::Info, true > LogVerbatim
constexpr std::array< uint8_t, layerIndexSize > layer
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:79
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 61 of file PatternRecognitionbyCLUE3D.cc.

References ticl::Trackster::charged_hadron, ticl::Trackster::electron, ticl::Trackster::neutral_hadron, pileupDistInMC::num, ticl::Trackster::photon, AlCaHLTBitMon_QueryRunRegistry::string, submitPVValidationJobs::t, and findQualityFiles::v.

63  {
65  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
66  << "[evt, tracksterId, cells, prob_photon, prob_ele, prob_chad, prob_nhad, layer_i, x_i, y_i, eta_i, phi_i, "
67  "energy_i, radius_i, rho_i, delta_i, isSeed_i";
68  }
69 
70  int num = 0;
71  const std::string sep(", ");
72  for (auto const &t : tracksters) {
73  for (auto v : t.vertices()) {
74  auto [lyrIdx, soaIdx] = layerIdx2layerandSoa[v];
75  auto const &thisLayer = clusters_[lyrIdx];
77  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
78  << "TracksterInfo: " << eventNumber << sep << num << sep << t.vertices().size() << sep
79  << t.id_probability(ticl::Trackster::ParticleType::photon) << sep
80  << t.id_probability(ticl::Trackster::ParticleType::electron) << sep
81  << t.id_probability(ticl::Trackster::ParticleType::charged_hadron) << sep
82  << t.id_probability(ticl::Trackster::ParticleType::neutral_hadron) << sep << lyrIdx << sep
83  << thisLayer.x[soaIdx] << sep << thisLayer.y[soaIdx] << sep << thisLayer.eta[soaIdx] << sep
84  << thisLayer.phi[soaIdx] << sep << thisLayer.energy[soaIdx] << sep << thisLayer.radius[soaIdx] << sep
85  << thisLayer.rho[soaIdx] << sep << thisLayer.delta[soaIdx] << sep << thisLayer.isSeed[soaIdx] << '\n';
86  }
87  }
88  num++;
89  }
90 }
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_
template<typename TILES >
void PatternRecognitionbyCLUE3D::energyRegressionAndID ( const std::vector< reco::CaloCluster > &  layerClusters,
const tensorflow::Session *  eidSession,
std::vector< Trackster > &  result 
)

Definition at line 289 of file PatternRecognitionbyCLUE3D.cc.

References a, funct::abs(), b, data, relval_parameters_module::energy, reco::CaloCluster::energy(), reco::CaloCluster::eta(), validate-o2o-wbm::f, lowptgsfeleseed::features(), reco::CaloCluster::hitsAndFractions(), mps_fire::i, input, dqmiolumiharvest::j, isotrackApplyRegressor::k, cmsLHEtoEOSManager::l, reco::CaloCluster::phi(), tensorflow::run(), and beam_dqm_sourceclient-live_cfg::vertices.

291  {
292  // Energy regression and particle identification strategy:
293  //
294  // 1. Set default values for regressed energy and particle id for each trackster.
295  // 2. Store indices of tracksters whose total sum of cluster energies is above the
296  // eidMinClusterEnergy_ (GeV) treshold. Inference is not applied for soft tracksters.
297  // 3. When no trackster passes the selection, return.
298  // 4. Create input and output tensors. The batch dimension is determined by the number of
299  // selected tracksters.
300  // 5. Fill input tensors with layer cluster features. Per layer, clusters are ordered descending
301  // by energy. Given that tensor data is contiguous in memory, we can use pointer arithmetic to
302  // fill values, even with batching.
303  // 6. Zero-fill features for empty clusters in each layer.
304  // 7. Batched inference.
305  // 8. Assign the regressed energy and id probabilities to each trackster.
306  //
307  // Indices used throughout this method:
308  // i -> batch element / trackster
309  // j -> layer
310  // k -> cluster
311  // l -> feature
312 
313  // set default values per trackster, determine if the cluster energy threshold is passed,
314  // and store indices of hard tracksters
315  std::vector<int> tracksterIndices;
316  for (int i = 0; i < static_cast<int>(tracksters.size()); i++) {
317  // calculate the cluster energy sum (2)
318  // note: after the loop, sumClusterEnergy might be just above the threshold which is enough to
319  // decide whether to run inference for the trackster or not
320  float sumClusterEnergy = 0.;
321  for (const unsigned int &vertex : tracksters[i].vertices()) {
322  sumClusterEnergy += static_cast<float>(layerClusters[vertex].energy());
323  // there might be many clusters, so try to stop early
324  if (sumClusterEnergy >= eidMinClusterEnergy_) {
325  // set default values (1)
326  tracksters[i].setRegressedEnergy(0.f);
327  tracksters[i].zeroProbabilities();
328  tracksterIndices.push_back(i);
329  break;
330  }
331  }
332  }
333 
334  // do nothing when no trackster passes the selection (3)
335  int batchSize = static_cast<int>(tracksterIndices.size());
336  if (batchSize == 0) {
337  return;
338  }
339 
340  // create input and output tensors (4)
341  tensorflow::TensorShape shape({batchSize, eidNLayers_, eidNClusters_, eidNFeatures_});
342  tensorflow::Tensor input(tensorflow::DT_FLOAT, shape);
344 
345  std::vector<tensorflow::Tensor> outputs;
346  std::vector<std::string> outputNames;
347  if (!eidOutputNameEnergy_.empty()) {
348  outputNames.push_back(eidOutputNameEnergy_);
349  }
350  if (!eidOutputNameId_.empty()) {
351  outputNames.push_back(eidOutputNameId_);
352  }
353 
354  // fill input tensor (5)
355  for (int i = 0; i < batchSize; i++) {
356  const Trackster &trackster = tracksters[tracksterIndices[i]];
357 
358  // per layer, we only consider the first eidNClusters_ clusters in terms of energy, so in order
359  // to avoid creating large / nested structures to do the sorting for an unknown number of total
360  // clusters, create a sorted list of layer cluster indices to keep track of the filled clusters
361  std::vector<int> clusterIndices(trackster.vertices().size());
362  for (int k = 0; k < (int)trackster.vertices().size(); k++) {
363  clusterIndices[k] = k;
364  }
365  sort(clusterIndices.begin(), clusterIndices.end(), [&layerClusters, &trackster](const int &a, const int &b) {
366  return layerClusters[trackster.vertices(a)].energy() > layerClusters[trackster.vertices(b)].energy();
367  });
368 
369  // keep track of the number of seen clusters per layer
370  std::vector<int> seenClusters(eidNLayers_);
371 
372  // loop through clusters by descending energy
373  for (const int &k : clusterIndices) {
374  // get features per layer and cluster and store the values directly in the input tensor
375  const reco::CaloCluster &cluster = layerClusters[trackster.vertices(k)];
376  int j = rhtools_.getLayerWithOffset(cluster.hitsAndFractions()[0].first) - 1;
377  if (j < eidNLayers_ && seenClusters[j] < eidNClusters_) {
378  // get the pointer to the first feature value for the current batch, layer and cluster
379  float *features = &input.tensor<float, 4>()(i, j, seenClusters[j], 0);
380 
381  // fill features
382  *(features++) = float(cluster.energy() / float(trackster.vertex_multiplicity(k)));
383  *(features++) = float(std::abs(cluster.eta()));
384  *(features) = float(cluster.phi());
385 
386  // increment seen clusters
387  seenClusters[j]++;
388  }
389  }
390 
391  // zero-fill features of empty clusters in each layer (6)
392  for (int j = 0; j < eidNLayers_; j++) {
393  for (int k = seenClusters[j]; k < eidNClusters_; k++) {
394  float *features = &input.tensor<float, 4>()(i, j, k, 0);
395  for (int l = 0; l < eidNFeatures_; l++) {
396  *(features++) = 0.f;
397  }
398  }
399  }
400  }
401 
402  // run the inference (7)
403  tensorflow::run(const_cast<tensorflow::Session *>(eidSession), inputList, outputNames, &outputs);
404 
405  // store regressed energy per trackster (8)
406  if (!eidOutputNameEnergy_.empty()) {
407  // get the pointer to the energy tensor, dimension is batch x 1
408  float *energy = outputs[0].flat<float>().data();
409 
410  for (const int &i : tracksterIndices) {
411  tracksters[i].setRegressedEnergy(*(energy++));
412  }
413  }
414 
415  // store id probabilities per trackster (8)
416  if (!eidOutputNameId_.empty()) {
417  // get the pointer to the id probability tensor, dimension is batch x id_probabilities.size()
418  int probsIdx = eidOutputNameEnergy_.empty() ? 0 : 1;
419  float *probs = outputs[probsIdx].flat<float>().data();
420 
421  for (const int &i : tracksterIndices) {
422  tracksters[i].setProbabilities(probs);
423  probs += tracksters[i].id_probabilities().size();
424  }
425  }
426 }
std::vector< NamedTensor > NamedTensorList
Definition: TensorFlow.h:30
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
static std::string const input
Definition: EdmProvDump.cc:47
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
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:213
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double energy() const
cluster energy
Definition: CaloCluster.h:149
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:365
double b
Definition: hdecay.h:118
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
double a
Definition: hdecay.h:119
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:184
template<typename TILES >
void PatternRecognitionbyCLUE3D::fillPSetDescription ( edm::ParameterSetDescription iDesc)
static

Definition at line 704 of file PatternRecognitionbyCLUE3D.cc.

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

704  {
705  iDesc.add<int>("algo_verbosity", 0);
706  iDesc.add<double>("criticalDensity", 4)->setComment("in GeV");
707  iDesc.add<int>("densitySiblingLayers", 3);
708  iDesc.add<double>("densityEtaPhiDistanceSqr", 0.0008);
709  iDesc.add<bool>("densityOnSameLayer", false);
710  iDesc.add<double>("criticalEtaPhiDistance", 0.035);
711  iDesc.add<double>("outlierMultiplier", 2);
712  iDesc.add<int>("minNumLayerCluster", 2)->setComment("Not Inclusive");
713  iDesc.add<std::string>("eid_input_name", "input");
714  iDesc.add<std::string>("eid_output_name_energy", "output/regressed_energy");
715  iDesc.add<std::string>("eid_output_name_id", "output/id_probabilities");
716  iDesc.add<double>("eid_min_cluster_energy", 1.);
717  iDesc.add<int>("eid_n_layers", 50);
718  iDesc.add<int>("eid_n_clusters", 10);
719 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
template<typename TILES >
int PatternRecognitionbyCLUE3D::findAndAssignTracksters ( const TILES &  tiles,
const std::vector< std::pair< int, int >> &  layerIdx2layerandSoa 
)
private

Definition at line 643 of file PatternRecognitionbyCLUE3D.cc.

References mps_fire::i, and phase1PixelTopology::layer.

644  {
645  unsigned int nTracksters = 0;
646 
647  std::vector<std::pair<int, int>> localStack;
648  // find cluster seeds and outlier
649  for (unsigned int layer = 0; layer < 2 * rhtools_.lastLayer(); layer++) {
650  auto &clustersOnLayer = clusters_[layer];
651  unsigned int numberOfClusters = clustersOnLayer.x.size();
652  for (unsigned int i = 0; i < numberOfClusters; i++) {
653  // initialize clusterIndex
654  clustersOnLayer.clusterIndex[i] = -1;
655  bool isSeed =
656  (clustersOnLayer.delta[i] > criticalEtaPhiDistance_) && (clustersOnLayer.rho[i] >= criticalDensity_);
657  bool isOutlier = (clustersOnLayer.delta[i] > outlierMultiplier_ * criticalEtaPhiDistance_) &&
658  (clustersOnLayer.rho[i] < criticalDensity_);
659  if (isSeed) {
661  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
662  << "Found seed on Layer " << layer << " SOAidx: " << i << " assigned ClusterIdx: " << nTracksters;
663  }
664  clustersOnLayer.clusterIndex[i] = nTracksters++;
665  clustersOnLayer.isSeed[i] = true;
666  localStack.emplace_back(layer, i);
667  } else if (!isOutlier) {
668  auto [lyrIdx, soaIdx] = clustersOnLayer.nearestHigher[i];
670  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
671  << "Found follower on Layer " << layer << " SOAidx: " << i << " attached to cluster on layer: " << lyrIdx
672  << " SOAidx: " << soaIdx;
673  }
674  if (lyrIdx >= 0)
675  clusters_[lyrIdx].followers[soaIdx].emplace_back(layer, i);
676  } else {
678  edm::LogVerbatim("PatternRecogntionbyCLUE3D")
679  << "Found Outlier on Layer " << layer << " SOAidx: " << i << " with rho: " << clustersOnLayer.rho[i]
680  << " and delta: " << clustersOnLayer.delta[i];
681  }
682  }
683  }
684  }
685 
686  // Propagate cluster index
687  while (!localStack.empty()) {
688  auto [lyrIdx, soaIdx] = localStack.back();
689  auto &thisSeed = clusters_[lyrIdx].followers[soaIdx];
690  localStack.pop_back();
691 
692  // loop over followers
693  for (auto [follower_lyrIdx, follower_soaIdx] : thisSeed) {
694  // pass id to a follower
695  clusters_[follower_lyrIdx].clusterIndex[follower_soaIdx] = clusters_[lyrIdx].clusterIndex[soaIdx];
696  // push this follower to localStack
697  localStack.emplace_back(follower_lyrIdx, follower_soaIdx);
698  }
699  }
700  return nTracksters;
701 }
Log< level::Info, true > LogVerbatim
std::vector< ClustersOnLayer > clusters_
constexpr std::array< uint8_t, layerIndexSize > layer
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:79
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 131 of file PatternRecognitionbyCLUE3D.cc.

References cms::cuda::assert(), ticl::assignPCAtoTracksters(), SplitLinear::begin, caloGeomToken_, dataset::end, ticl::PatternRecognitionAlgoBaseT< TILES >::Inputs::es, ticl::PatternRecognitionAlgoBaseT< TILES >::Inputs::ev, edm::EventAuxiliary::event(), edm::Event::eventAuxiliary(), validate-o2o-wbm::f, relativeConstraints::geom, edm::EventSetup::getData(), mps_fire::i, phase1PixelTopology::layer, ticl::PatternRecognitionAlgoBaseT< TILES >::Inputs::layerClusters, ticl::PatternRecognitionAlgoBaseT< TILES >::Inputs::layerClustersTime, ticl::PatternRecognitionAlgoBaseT< TILES >::Inputs::mask, SiStripPI::max, point, ticl::PatternRecognitionAlgoBaseT< TILES >::Inputs::regions, reset(), mathSSE::sqrt(), makeMuonMisalignmentScenario::sum_x, makeMuonMisalignmentScenario::sum_y, submitPVValidationJobs::t, ticl::PatternRecognitionAlgoBaseT< TILES >::Inputs::tfSession, ticl::PatternRecognitionAlgoBaseT< TILES >::Inputs::tiles, and findQualityFiles::v.

134  {
135  // Protect from events with no seeding regions
136  if (input.regions.empty())
137  return;
138 
139  const int eventNumber = input.ev.eventAuxiliary().event();
141  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "New Event";
142  }
143 
144  edm::EventSetup const &es = input.es;
146  rhtools_.setGeometry(geom);
147 
148  clusters_.clear();
149  clusters_.resize(2 * rhtools_.lastLayer(false));
150  std::vector<std::pair<int, int>> layerIdx2layerandSoa; //used everywhere also to propagate cluster masking
151 
152  layerIdx2layerandSoa.reserve(input.layerClusters.size());
153  unsigned int layerIdx = 0;
154  for (auto const &lc : input.layerClusters) {
155  if (input.mask[layerIdx] == 0.) {
157  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Skipping masked clustrer: " << layerIdx;
158  }
159  layerIdx2layerandSoa.emplace_back(-1, -1);
160  layerIdx++;
161  continue;
162  }
163  const auto firstHitDetId = lc.hitsAndFractions()[0].first;
164  int layer = rhtools_.getLayerWithOffset(firstHitDetId) - 1 +
165  rhtools_.lastLayer(false) * ((rhtools_.zside(firstHitDetId) + 1) >> 1);
166  assert(layer >= 0);
167 
168  layerIdx2layerandSoa.emplace_back(layer, clusters_[layer].x.size());
169  float sum_x = 0.;
170  float sum_y = 0.;
171  float sum_sqr_x = 0.;
172  float sum_sqr_y = 0.;
173  float ref_x = lc.x();
174  float ref_y = lc.y();
175  float invClsize = 1. / lc.hitsAndFractions().size();
176  for (auto const &hitsAndFractions : lc.hitsAndFractions()) {
177  auto const &point = rhtools_.getPosition(hitsAndFractions.first);
178  sum_x += point.x() - ref_x;
179  sum_sqr_x += (point.x() - ref_x) * (point.x() - ref_x);
180  sum_y += point.y() - ref_y;
181  sum_sqr_y += (point.y() - ref_y) * (point.y() - ref_y);
182  }
183  // The variance of X for X uniform in circle of radius R is R^2/4,
184  // therefore we multiply the sqrt(var) by 2 to have a rough estimate of the
185  // radius. On the other hand, while averaging the x and y radius, we would
186  // end up dividing by 2. Hence we omit the value here and in the average
187  // below, too.
188  float radius_x = sqrt((sum_sqr_x - (sum_x * sum_x) * invClsize) * invClsize);
189  float radius_y = sqrt((sum_sqr_y - (sum_y * sum_y) * invClsize) * invClsize);
190  clusters_[layer].x.emplace_back(lc.x());
191  clusters_[layer].y.emplace_back(lc.y());
192  clusters_[layer].radius.emplace_back(radius_x + radius_y);
193  clusters_[layer].eta.emplace_back(lc.eta());
194  clusters_[layer].phi.emplace_back(lc.phi());
195  clusters_[layer].cells.push_back(lc.hitsAndFractions().size());
196  clusters_[layer].energy.emplace_back(lc.energy());
197  clusters_[layer].isSeed.push_back(false);
198  clusters_[layer].clusterIndex.emplace_back(-1);
199  clusters_[layer].layerClusterOriginalIdx.emplace_back(layerIdx++);
200  clusters_[layer].nearestHigher.emplace_back(-1, -1);
201  clusters_[layer].rho.emplace_back(0.f);
202  clusters_[layer].delta.emplace_back(std::numeric_limits<float>::max());
203  }
204  for (unsigned int layer = 0; layer < clusters_.size(); layer++) {
205  clusters_[layer].followers.resize(clusters_[layer].x.size());
206  }
207 
208  auto lastLayerPerSide = static_cast<int>(rhtools_.lastLayer(false));
209  int maxLayer = 2 * lastLayerPerSide - 1;
210  std::vector<int> numberOfClustersPerLayer(maxLayer, 0);
211  for (int i = 0; i <= maxLayer; i++) {
212  calculateLocalDensity(input.tiles, i, layerIdx2layerandSoa);
213  }
214  for (int i = 0; i <= maxLayer; i++) {
215  calculateDistanceToHigher(input.tiles, i, layerIdx2layerandSoa);
216  }
217 
218  auto nTracksters = findAndAssignTracksters(input.tiles, layerIdx2layerandSoa);
220  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Reconstructed " << nTracksters << " tracksters" << std::endl;
221  dumpClusters(layerIdx2layerandSoa, eventNumber);
222  }
223 
224  // Build Trackster
225  result.resize(nTracksters);
226 
227  for (unsigned int layer = 0; layer < clusters_.size(); ++layer) {
228  const auto &thisLayer = clusters_[layer];
230  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Examining Layer: " << layer;
231  }
232  for (unsigned int lc = 0; lc < thisLayer.x.size(); ++lc) {
234  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Trackster " << thisLayer.clusterIndex[lc];
235  }
236  if (thisLayer.clusterIndex[lc] >= 0) {
238  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << " adding lcIdx: " << thisLayer.layerClusterOriginalIdx[lc];
239  }
240  result[thisLayer.clusterIndex[lc]].vertices().push_back(thisLayer.layerClusterOriginalIdx[lc]);
241  result[thisLayer.clusterIndex[lc]].vertex_multiplicity().push_back(1);
242  // loop over followers
243  for (auto [follower_lyrIdx, follower_soaIdx] : thisLayer.followers[lc]) {
244  std::array<unsigned int, 2> edge = {
245  {(unsigned int)thisLayer.layerClusterOriginalIdx[lc],
246  (unsigned int)clusters_[follower_lyrIdx].layerClusterOriginalIdx[follower_soaIdx]}};
247  result[thisLayer.clusterIndex[lc]].edges().push_back(edge);
248  }
249  }
250  }
251  }
252 
253  result.erase(
254  std::remove_if(std::begin(result),
255  std::end(result),
256  [&](auto const &v) { return static_cast<int>(v.vertices().size()) < minNumLayerCluster_; }),
257  result.end());
258  result.shrink_to_fit();
259 
261  input.layerClusters,
262  input.layerClustersTime,
263  rhtools_.getPositionLayer(rhtools_.lastLayerEE(false), false).z());
264 
265  // run energy regression and ID
266  energyRegressionAndID(input.layerClusters, input.tfSession, result);
268  for (auto const &t : result) {
269  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Barycenter: " << t.barycenter();
270  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "LCs: " << t.vertices().size();
271  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Energy: " << t.raw_energy();
272  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << "Regressed: " << t.regressed_energy();
273  }
274  }
275 
276  // Dump Tracksters information
278  dumpTracksters(layerIdx2layerandSoa, eventNumber, result);
279  }
280 
281  // Reset internal clusters_ structure of array for next event
282  reset();
284  edm::LogVerbatim("PatternRecogntionbyCLUE3D") << std::endl;
285  }
286 }
Log< level::Info, true > LogVerbatim
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
void dumpTracksters(const std::vector< std::pair< int, int >> &layerIdx2layerandSoa, const int, const std::vector< Trackster > &) const
int zside(const DetId &id) const
Definition: RecHitTools.cc:163
void energyRegressionAndID(const std::vector< reco::CaloCluster > &layerClusters, const tensorflow::Session *, std::vector< Trackster > &result)
std::vector< ClustersOnLayer > clusters_
assert(be >=bs)
void calculateDistanceToHigher(const TILES &, const int layerId, const std::vector< std::pair< int, int >> &)
constexpr std::array< uint8_t, layerIndexSize > layer
static std::string const input
Definition: EdmProvDump.cc:47
tuple result
Definition: mps_fire.py:311
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void assignPCAtoTracksters(std::vector< Trackster > &, const std::vector< reco::CaloCluster > &, const edm::ValueMap< std::pair< float, float >> &, double, bool energyWeight=true)
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:79
int findAndAssignTracksters(const TILES &, const std::vector< std::pair< int, int >> &)
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:365
GlobalPoint getPositionLayer(int layer, bool nose=false) const
Definition: RecHitTools.cc:141
void dumpClusters(const std::vector< std::pair< int, int >> &layerIdx2layerandSoa, const int) const
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:68
void calculateLocalDensity(const TILES &, const int layerId, const std::vector< std::pair< int, int >> &)
unsigned int lastLayerEE(bool nose=false) const
Definition: RecHitTools.h:75
string end
Definition: dataset.py:937
*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
template<typename TILES >
void ticl::PatternRecognitionbyCLUE3D< TILES >::reset ( void  )
inlineprivate

Definition at line 81 of file PatternRecognitionbyCLUE3D.h.

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

81  {
82  for (auto& c : clusters_) {
83  c.clear();
84  c.shrink_to_fit();
85  }
86  }
const edm::EventSetup & c
std::vector< ClustersOnLayer > clusters_

Member Data Documentation

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

Definition at line 98 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 99 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 103 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 101 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 102 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 100 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 107 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 110 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 112 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 117 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 111 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 108 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 109 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 115 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 106 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 105 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 104 of file PatternRecognitionbyCLUE3D.h.

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

Definition at line 114 of file PatternRecognitionbyCLUE3D.h.