7 #include "tbb/task_arena.h"
23 template <
typename TILES>
29 criticalDensity_(conf.getParameter<double>(
"criticalDensity")),
30 densitySiblingLayers_(conf.getParameter<int>(
"densitySiblingLayers")),
31 densityEtaPhiDistanceSqr_(conf.getParameter<double>(
"densityEtaPhiDistanceSqr")),
32 densityOnSameLayer_(conf.getParameter<bool>(
"densityOnSameLayer")),
33 criticalEtaPhiDistance_(conf.getParameter<double>(
"criticalEtaPhiDistance")),
34 outlierMultiplier_(conf.getParameter<double>(
"outlierMultiplier")),
35 minNumLayerCluster_(conf.getParameter<int>(
"minNumLayerCluster")),
36 eidInputName_(conf.getParameter<std::
string>(
"eid_input_name")),
37 eidOutputNameEnergy_(conf.getParameter<std::
string>(
"eid_output_name_energy")),
38 eidOutputNameId_(conf.getParameter<std::
string>(
"eid_output_name_id")),
39 eidMinClusterEnergy_(conf.getParameter<double>(
"eid_min_cluster_energy")),
40 eidNLayers_(conf.getParameter<int>(
"eid_n_layers")),
41 eidNClusters_(conf.getParameter<int>(
"eid_n_clusters")),
42 eidSession_(nullptr) {
45 if (trackstersCache ==
nullptr || trackstersCache->
eidGraphDef ==
nullptr) {
47 <<
"PatternRecognitionbyCLUE3D received an empty graph definition from the global cache";
52 template <
typename TILES>
56 auto lastLayerPerSide = (
unsigned int)(rhtools_.lastLayer(
false));
57 unsigned int maxLayer = 2 * lastLayerPerSide - 1;
59 for (
int ieta = 0; ieta < nEtaBin; ieta++) {
60 auto offset = ieta * nPhiBin;
61 for (
int phi = 0; phi < nPhiBin; phi++) {
62 int iphi = ((phi % nPhiBin + nPhiBin) % nPhiBin);
64 if (this->algo_verbosity_ > this->Advanced) {
65 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"Layer: " <<
layer <<
" ieta: " << ieta <<
" phi: " << phi
74 template <
typename TILES>
76 const int eventNumber,
77 const std::vector<Trackster> &tracksters)
const {
80 <<
"[evt, tracksterId, cells, prob_photon, prob_ele, prob_chad, prob_nhad, layer_i, x_i, y_i, eta_i, phi_i, "
81 "energy_i, radius_i, rho_i, delta_i, isSeed_i";
86 for (
auto const &
t : tracksters) {
87 for (
auto v :
t.vertices()) {
88 auto [lyrIdx, soaIdx] = layerIdx2layerandSoa[
v];
89 auto const &thisLayer = clusters_[lyrIdx];
92 <<
"TracksterInfo: " << eventNumber << sep << num << sep <<
t.vertices().size() << sep
97 << thisLayer.x[soaIdx] << sep << thisLayer.y[soaIdx] << sep << thisLayer.eta[soaIdx] << sep
98 << thisLayer.phi[soaIdx] << sep << thisLayer.energy[soaIdx] << sep << thisLayer.radius[soaIdx] << sep
99 << thisLayer.rho[soaIdx] << sep << thisLayer.delta[soaIdx] << sep << thisLayer.isSeed[soaIdx] <<
'\n';
106 template <
typename TILES>
108 const int eventNumber)
const {
110 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"[evt, layer, x, y, eta, phi, cells, energy, radius, rho, delta, "
111 "isSeed, clusterIdx, layerClusterOriginalIdx";
118 auto const &thisLayer = clusters_[
layer];
120 for (
auto v : thisLayer.x) {
123 <<
"ClusterInfo: " << eventNumber <<
", " <<
layer <<
", " <<
v <<
", " << thisLayer.y[
num] <<
", "
124 << thisLayer.eta[
num] <<
", " << thisLayer.phi[
num] <<
", " << thisLayer.cells[
num] <<
", "
125 << thisLayer.energy[
num] <<
", " << thisLayer.radius[
num] <<
", " << thisLayer.rho[
num] <<
", "
126 << thisLayer.delta[
num] <<
", " << thisLayer.isSeed[
num] <<
", " << thisLayer.clusterIndex[
num] <<
", "
127 << thisLayer.layerClusterOriginalIdx[
num];
132 for (
unsigned int lcIdx = 0; lcIdx < layerIdx2layerandSoa.size(); lcIdx++) {
133 auto const &layerandSoa = layerIdx2layerandSoa[lcIdx];
135 if ((layerandSoa.first == -1) && (layerandSoa.second == -1))
139 <<
"lcIdx: " << lcIdx <<
" on Layer: " << layerandSoa.first <<
" SOA: " << layerandSoa.second;
144 template <
typename TILES>
147 std::vector<Trackster> &
result,
148 std::unordered_map<
int, std::vector<int>> &seedToTracksterAssociation) {
160 rhtools_.setGeometry(geom);
163 clusters_.resize(2 * rhtools_.lastLayer(
false));
164 std::vector<std::pair<int, int>> layerIdx2layerandSoa;
167 unsigned int layerIdx = 0;
169 if (input.
mask[layerIdx] == 0.) {
171 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"Skipping masked clustrer: " << layerIdx;
173 layerIdx2layerandSoa.emplace_back(-1, -1);
177 const auto firstHitDetId = lc.hitsAndFractions()[0].first;
178 int layer = rhtools_.getLayerWithOffset(firstHitDetId) - 1 +
179 rhtools_.lastLayer(
false) * ((rhtools_.zside(firstHitDetId) + 1) >> 1);
182 layerIdx2layerandSoa.emplace_back(layer, clusters_[layer].x.size());
185 float sum_sqr_x = 0.;
186 float sum_sqr_y = 0.;
187 float ref_x = lc.x();
188 float ref_y = lc.y();
189 float invClsize = 1. / lc.hitsAndFractions().size();
190 for (
auto const &hitsAndFractions : lc.hitsAndFractions()) {
191 auto const &
point = rhtools_.getPosition(hitsAndFractions.first);
192 sum_x +=
point.x() - ref_x;
193 sum_sqr_x += (
point.x() - ref_x) * (
point.x() - ref_x);
194 sum_y +=
point.y() - ref_y;
195 sum_sqr_y += (
point.y() - ref_y) * (
point.y() - ref_y);
202 float radius_x =
sqrt((sum_sqr_x - (sum_x * sum_x) * invClsize) * invClsize);
203 float radius_y =
sqrt((sum_sqr_y - (sum_y * sum_y) * invClsize) * invClsize);
204 clusters_[
layer].x.emplace_back(lc.x());
205 clusters_[
layer].y.emplace_back(lc.y());
206 clusters_[
layer].radius.emplace_back(radius_x + radius_y);
207 clusters_[
layer].eta.emplace_back(lc.eta());
208 clusters_[
layer].phi.emplace_back(lc.phi());
209 clusters_[
layer].cells.push_back(lc.hitsAndFractions().size());
210 clusters_[
layer].energy.emplace_back(lc.energy());
211 clusters_[
layer].isSeed.push_back(
false);
212 clusters_[
layer].clusterIndex.emplace_back(-1);
213 clusters_[
layer].layerClusterOriginalIdx.emplace_back(layerIdx++);
214 clusters_[
layer].nearestHigher.emplace_back(-1, -1);
215 clusters_[
layer].rho.emplace_back(0.
f);
219 clusters_[
layer].followers.resize(clusters_[
layer].x.size());
222 auto lastLayerPerSide = (
unsigned int)(rhtools_.lastLayer(
false));
223 unsigned int maxLayer = 2 * lastLayerPerSide - 1;
224 std::vector<int> numberOfClustersPerLayer(maxLayer, 0);
225 for (
unsigned int i = 0;
i <= maxLayer;
i++) {
226 calculateLocalDensity(input.
tiles,
i, layerIdx2layerandSoa);
228 for (
unsigned int i = 0;
i <= maxLayer;
i++) {
229 calculateDistanceToHigher(input.
tiles,
i, layerIdx2layerandSoa);
232 auto nTracksters = findAndAssignTracksters(input.
tiles, layerIdx2layerandSoa);
234 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"Reconstructed " << nTracksters <<
" tracksters" << std::endl;
235 dumpClusters(layerIdx2layerandSoa, eventNumber);
239 result.resize(nTracksters);
242 const auto &thisLayer = clusters_[
layer];
246 for (
unsigned int lc = 0; lc < thisLayer.x.size(); ++lc) {
248 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"Trackster " << thisLayer.clusterIndex[lc];
250 if (thisLayer.clusterIndex[lc] >= 0) {
252 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
" adding lcIdx: " << thisLayer.layerClusterOriginalIdx[lc];
254 result[thisLayer.clusterIndex[lc]].vertices().push_back(thisLayer.layerClusterOriginalIdx[lc]);
255 result[thisLayer.clusterIndex[lc]].vertex_multiplicity().push_back(1);
257 for (
auto [follower_lyrIdx, follower_soaIdx] : thisLayer.followers[lc]) {
258 std::array<unsigned int, 2> edge = {
259 {(
unsigned int)thisLayer.layerClusterOriginalIdx[lc],
260 (
unsigned int)clusters_[follower_lyrIdx].layerClusterOriginalIdx[follower_soaIdx]}};
261 result[thisLayer.clusterIndex[lc]].edges().push_back(edge);
270 [&](
auto const &
v) {
return static_cast<int>(
v.vertices().size()) < minNumLayerCluster_; }),
272 result.shrink_to_fit();
277 rhtools_.getPositionLayer(rhtools_.lastLayerEE(
false),
false).z());
282 for (
auto const &
t : result) {
283 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"Barycenter: " <<
t.barycenter();
286 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"Regressed: " <<
t.regressed_energy();
292 dumpTracksters(layerIdx2layerandSoa, eventNumber, result);
302 template <
typename TILES>
304 std::vector<Trackster> &tracksters) {
328 std::vector<int> tracksterIndices;
329 for (
int i = 0; i < static_cast<int>(tracksters.size());
i++) {
333 float sumClusterEnergy = 0.;
334 for (
const unsigned int &vertex : tracksters[
i].
vertices()) {
335 sumClusterEnergy +=
static_cast<float>(layerClusters[vertex].energy());
337 if (sumClusterEnergy >= eidMinClusterEnergy_) {
339 tracksters[
i].setRegressedEnergy(0.
f);
340 tracksters[
i].zeroProbabilities();
341 tracksterIndices.push_back(
i);
348 int batchSize =
static_cast<int>(tracksterIndices.size());
349 if (batchSize == 0) {
354 tensorflow::TensorShape shape({batchSize, eidNLayers_, eidNClusters_, eidNFeatures_});
355 tensorflow::Tensor
input(tensorflow::DT_FLOAT, shape);
358 std::vector<tensorflow::Tensor> outputs;
359 std::vector<std::string> outputNames;
360 if (!eidOutputNameEnergy_.empty()) {
361 outputNames.push_back(eidOutputNameEnergy_);
363 if (!eidOutputNameId_.empty()) {
364 outputNames.push_back(eidOutputNameId_);
368 for (
int i = 0;
i < batchSize;
i++) {
369 const Trackster &trackster = tracksters[tracksterIndices[
i]];
374 std::vector<int> clusterIndices(trackster.vertices().size());
375 for (
int k = 0;
k < (int)trackster.vertices().size();
k++) {
376 clusterIndices[
k] =
k;
378 sort(clusterIndices.begin(), clusterIndices.end(), [&layerClusters, &trackster](
const int &
a,
const int &
b) {
379 return layerClusters[trackster.vertices(a)].energy() > layerClusters[trackster.vertices(
b)].energy();
383 std::vector<int> seenClusters(eidNLayers_);
386 for (
const int &
k : clusterIndices) {
390 if (j < eidNLayers_ && seenClusters[j] < eidNClusters_) {
392 float *
features = &input.tensor<float, 4>()(
i, j, seenClusters[j], 0);
395 *(features++) =
float(cluster.
energy() / float(trackster.vertex_multiplicity(
k)));
405 for (
int j = 0;
j < eidNLayers_;
j++) {
406 for (
int k = seenClusters[
j];
k < eidNClusters_;
k++) {
407 float *
features = &input.tensor<float, 4>()(
i,
j,
k, 0);
408 for (
int l = 0;
l < eidNFeatures_;
l++) {
419 if (!eidOutputNameEnergy_.empty()) {
421 float *
energy = outputs[0].flat<
float>().
data();
423 for (
const int &
i : tracksterIndices) {
424 tracksters[
i].setRegressedEnergy(*(energy++));
429 if (!eidOutputNameId_.empty()) {
431 int probsIdx = eidOutputNameEnergy_.empty() ? 0 : 1;
432 float *probs = outputs[probsIdx].flat<
float>().
data();
434 for (
const int &
i : tracksterIndices) {
435 tracksters[
i].setProbabilities(probs);
436 probs += tracksters[
i].id_probabilities().size();
441 template <
typename TILES>
443 const TILES &tiles,
const unsigned int layerId,
const std::vector<std::pair<int, int>> &layerIdx2layerandSoa) {
446 auto &clustersOnLayer = clusters_[layerId];
447 unsigned int numberOfClusters = clustersOnLayer.x.size();
449 auto isReachable = [&](
float x1,
float x2,
float y1,
float y2,
float delta_sqr) ->
bool {
452 << ((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2)) <<
" vs " << delta_sqr <<
"["
453 << ((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2) < delta_sqr) <<
"]\n";
455 return (x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2) < delta_sqr;
458 for (
unsigned int i = 0;
i < numberOfClusters;
i++) {
460 auto lastLayerPerSide =
static_cast<unsigned int>(rhtools_.lastLayer(
false)) - 1;
462 unsigned int maxLayer = 2 * lastLayerPerSide - 1;
463 if (layerId < lastLayerPerSide) {
464 minLayer =
std::max(layerId - densitySiblingLayers_, minLayer);
465 maxLayer =
std::min(layerId + densitySiblingLayers_, lastLayerPerSide - 1);
467 minLayer =
std::max(layerId - densitySiblingLayers_, lastLayerPerSide + 1);
468 maxLayer =
std::min(layerId + densitySiblingLayers_, maxLayer);
470 for (
unsigned int currentLayer = minLayer; currentLayer <= maxLayer; currentLayer++) {
472 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"RefLayer: " << layerId <<
" SoaIDX: " <<
i;
473 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"NextLayer: " << currentLayer;
475 const auto &tileOnLayer = tiles[currentLayer];
476 bool onSameLayer = (currentLayer == layerId);
478 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"onSameLayer: " << onSameLayer;
480 int etaWindow = onSameLayer ? 2 : 1;
481 int phiWindow = onSameLayer ? 2 : 1;
482 int etaBinMin =
std::max(tileOnLayer.etaBin(clustersOnLayer.eta[
i]) - etaWindow, 0);
483 int etaBinMax =
std::min(tileOnLayer.etaBin(clustersOnLayer.eta[
i]) + etaWindow, nEtaBin);
484 int phiBinMin = tileOnLayer.phiBin(clustersOnLayer.phi[
i]) - phiWindow;
485 int phiBinMax = tileOnLayer.phiBin(clustersOnLayer.phi[
i]) + phiWindow;
487 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"eta: " << clustersOnLayer.eta[
i];
488 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"phi: " << clustersOnLayer.phi[
i];
489 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"etaBinMin: " << etaBinMin <<
", etaBinMax: " << etaBinMax;
490 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"phiBinMin: " << phiBinMin <<
", phiBinMax: " << phiBinMax;
492 for (
int ieta = etaBinMin; ieta <= etaBinMax; ++ieta) {
493 auto offset = ieta * nPhiBin;
497 for (
int iphi_it = phiBinMin; iphi_it <= phiBinMax; ++iphi_it) {
498 int iphi = ((iphi_it % nPhiBin + nPhiBin) % nPhiBin);
502 <<
"Entries in tileBin: " << tileOnLayer[
offset + iphi].size();
504 for (
auto otherClusterIdx : tileOnLayer[
offset + iphi]) {
505 auto const &layerandSoa = layerIdx2layerandSoa[otherClusterIdx];
507 if ((layerandSoa.first == -1) && (layerandSoa.second == -1)) {
509 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"Skipping masked layerIdx " << otherClusterIdx;
513 auto const &clustersLayer = clusters_[layerandSoa.first];
516 <<
"OtherLayer: " << layerandSoa.first <<
" SoaIDX: " << layerandSoa.second;
517 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"OtherEta: " << clustersLayer.eta[layerandSoa.second];
518 edm::LogVerbatim(
"PatternRecogntionbyCLUE3D") <<
"OtherPhi: " << clustersLayer.phi[layerandSoa.second];
521 float delta = clustersOnLayer.radius[
i] + clustersLayer.radius[layerandSoa.second] + 2.6f;
522 if (densityOnSameLayer_ && onSameLayer) {
523 if (isReachable(clustersOnLayer.x[
i],
524 clustersLayer.x[layerandSoa.second],
525 clustersOnLayer.y[
i],
526 clustersLayer.y[layerandSoa.second],
528 clustersOnLayer.rho[
i] += (clustersOnLayer.layerClusterOriginalIdx[
i] == otherClusterIdx ? 1.f : 0.2f) *
529 clustersLayer.energy[layerandSoa.second];
535 clustersOnLayer.phi[
i],
536 clustersLayer.eta[layerandSoa.second],
537 clustersLayer.phi[layerandSoa.second]);
540 clustersOnLayer.phi[
i],
541 clustersLayer.eta[layerandSoa.second],
542 clustersLayer.phi[layerandSoa.second]) < densityEtaPhiDistanceSqr_) {
543 auto energyToAdd = (clustersOnLayer.layerClusterOriginalIdx[
i] == otherClusterIdx ? 1.f : 0.5f) *
544 clustersLayer.energy[layerandSoa.second];
545 clustersOnLayer.rho[
i] += energyToAdd;
547 <<
"Adding " << energyToAdd <<
" partial " << clustersOnLayer.rho[
i];
560 template <
typename TILES>
562 const TILES &tiles,
const unsigned int layerId,
const std::vector<std::pair<int, int>> &layerIdx2layerandSoa) {
565 auto &clustersOnLayer = clusters_[layerId];
566 unsigned int numberOfClusters = clustersOnLayer.x.size();
568 for (
unsigned int i = 0;
i < numberOfClusters;
i++) {
571 <<
"Starting searching nearestHigher on " << layerId <<
" with rho: " << clustersOnLayer.rho[
i]
572 <<
" at eta, phi: " << tiles[layerId].etaBin(clustersOnLayer.eta[
i]) <<
", "
573 << tiles[layerId].etaBin(clustersOnLayer.phi[
i]);
576 auto lastLayerPerSide =
static_cast<unsigned int>(rhtools_.lastLayer(
false));
578 unsigned int maxLayer = 2 * lastLayerPerSide - 1;
579 if (layerId < lastLayerPerSide) {
580 minLayer =
std::max(layerId - densitySiblingLayers_, minLayer);
581 maxLayer =
std::min(layerId + densitySiblingLayers_, lastLayerPerSide - 1);
583 minLayer =
std::max(layerId - densitySiblingLayers_, lastLayerPerSide + 1);
584 maxLayer =
std::min(layerId + densitySiblingLayers_, maxLayer);
587 float i_delta = maxDelta;
588 std::pair<int, int> i_nearestHigher(-1, -1);
589 for (
unsigned int currentLayer = minLayer; currentLayer <= maxLayer; currentLayer++) {
590 const auto &tileOnLayer = tiles[currentLayer];
593 int etaBinMin =
std::max(tileOnLayer.etaBin(clustersOnLayer.eta[
i]) - etaWindow, 0);
594 int etaBinMax =
std::min(tileOnLayer.etaBin(clustersOnLayer.eta[
i]) + etaWindow, nEtaBin);
595 int phiBinMin = tileOnLayer.phiBin(clustersOnLayer.phi[
i]) - phiWindow;
596 int phiBinMax = tileOnLayer.phiBin(clustersOnLayer.phi[
i]) + phiWindow;
597 for (
int ieta = etaBinMin; ieta <= etaBinMax; ++ieta) {
598 auto offset = ieta * nPhiBin;
599 for (
int iphi_it = phiBinMin; iphi_it <= phiBinMax; ++iphi_it) {
600 int iphi = ((iphi_it % nPhiBin + nPhiBin) % nPhiBin);
603 <<
"Searching nearestHigher on " << currentLayer <<
" eta, phi: " << ieta <<
", " << iphi_it;
605 for (
auto otherClusterIdx : tileOnLayer[
offset + iphi]) {
606 auto const &layerandSoa = layerIdx2layerandSoa[otherClusterIdx];
608 if ((layerandSoa.first == -1) && (layerandSoa.second == -1))
610 auto const &clustersOnOtherLayer = clusters_[layerandSoa.first];
612 clustersOnLayer.phi[
i],
613 clustersOnOtherLayer.eta[layerandSoa.second],
614 clustersOnOtherLayer.phi[layerandSoa.second]);
615 bool foundHigher = (clustersOnOtherLayer.rho[layerandSoa.second] > clustersOnLayer.rho[
i]) ||
616 (clustersOnOtherLayer.rho[layerandSoa.second] == clustersOnLayer.rho[
i] &&
617 clustersOnOtherLayer.layerClusterOriginalIdx[layerandSoa.second] >
618 clustersOnLayer.layerClusterOriginalIdx[
i]);
621 <<
"Searching nearestHigher on " << currentLayer
622 <<
" with rho: " << clustersOnOtherLayer.rho[layerandSoa.second]
623 <<
" on layerIdxInSOA: " << layerandSoa.first <<
", " << layerandSoa.second
624 <<
" with distance: " <<
sqrt(dist) <<
" foundHigher: " << foundHigher;
626 if (foundHigher && dist <= i_delta) {
630 i_nearestHigher = layerandSoa;
637 bool foundNearestHigherInEtaPhiCylinder = (i_delta != maxDelta);
640 <<
"i_delta: " <<
sqrt(i_delta) <<
" passed: " << foundNearestHigherInEtaPhiCylinder <<
" "
641 << i_nearestHigher.first <<
" " << i_nearestHigher.second;
643 if (foundNearestHigherInEtaPhiCylinder) {
644 clustersOnLayer.delta[
i] =
sqrt(i_delta);
645 clustersOnLayer.nearestHigher[
i] = i_nearestHigher;
649 clustersOnLayer.delta[
i] = maxDelta;
650 clustersOnLayer.nearestHigher[
i] = {-1, -1};
655 template <
typename TILES>
657 const TILES &tiles,
const std::vector<std::pair<int, int>> &layerIdx2layerandSoa) {
658 unsigned int nTracksters = 0;
660 std::vector<std::pair<int, int>> localStack;
662 for (
unsigned int layer = 0;
layer < 2 * rhtools_.lastLayer();
layer++) {
663 auto &clustersOnLayer = clusters_[
layer];
664 unsigned int numberOfClusters = clustersOnLayer.x.size();
665 for (
unsigned int i = 0;
i < numberOfClusters;
i++) {
667 clustersOnLayer.clusterIndex[
i] = -1;
669 (clustersOnLayer.delta[
i] > criticalEtaPhiDistance_) && (clustersOnLayer.rho[
i] >= criticalDensity_);
670 bool isOutlier = (clustersOnLayer.delta[
i] > outlierMultiplier_ * criticalEtaPhiDistance_) &&
671 (clustersOnLayer.rho[
i] < criticalDensity_);
675 <<
"Found seed on Layer " <<
layer <<
" SOAidx: " <<
i <<
" assigned ClusterIdx: " << nTracksters;
677 clustersOnLayer.clusterIndex[
i] = nTracksters++;
678 clustersOnLayer.isSeed[
i] =
true;
679 localStack.emplace_back(
layer,
i);
680 }
else if (!isOutlier) {
681 auto [lyrIdx, soaIdx] = clustersOnLayer.nearestHigher[
i];
684 <<
"Found follower on Layer " <<
layer <<
" SOAidx: " <<
i <<
" attached to cluster on layer: " << lyrIdx
685 <<
" SOAidx: " << soaIdx;
687 clusters_[lyrIdx].followers[soaIdx].emplace_back(
layer,
i);
691 <<
"Found Outlier on Layer " <<
layer <<
" SOAidx: " <<
i <<
" with rho: " << clustersOnLayer.rho[
i]
692 <<
" and delta: " << clustersOnLayer.delta[
i];
699 while (!localStack.empty()) {
700 auto [lyrIdx, soaIdx] = localStack.back();
701 auto &thisSeed = clusters_[lyrIdx].followers[soaIdx];
702 localStack.pop_back();
705 for (
auto [follower_lyrIdx, follower_soaIdx] : thisSeed) {
707 clusters_[follower_lyrIdx].clusterIndex[follower_soaIdx] = clusters_[lyrIdx].clusterIndex[soaIdx];
709 localStack.emplace_back(follower_lyrIdx, follower_soaIdx);
715 template <
typename TILES>
717 iDesc.
add<
int>(
"algo_verbosity", 0);
718 iDesc.
add<
double>(
"criticalDensity", 4)->setComment(
"in GeV");
719 iDesc.
add<
int>(
"densitySiblingLayers", 3);
720 iDesc.
add<
double>(
"densityEtaPhiDistanceSqr", 0.0008);
721 iDesc.
add<
bool>(
"densityOnSameLayer",
false);
722 iDesc.
add<
double>(
"criticalEtaPhiDistance", 0.035);
723 iDesc.
add<
double>(
"outlierMultiplier", 2);
724 iDesc.
add<
int>(
"minNumLayerCluster", 5)->setComment(
"Not Inclusive");
726 iDesc.
add<
std::string>(
"eid_output_name_energy",
"output/regressed_energy");
727 iDesc.
add<
std::string>(
"eid_output_name_id",
"output/id_probabilities");
728 iDesc.
add<
double>(
"eid_min_cluster_energy", 1.);
729 iDesc.
add<
int>(
"eid_n_layers", 50);
730 iDesc.
add<
int>(
"eid_n_clusters", 10);
Session * createSession(SessionOptions &sessionOptions)
Log< level::Info, true > LogVerbatim
std::vector< NamedTensor > NamedTensorList
void dumpTracksters(const std::vector< std::pair< int, int >> &layerIdx2layerandSoa, const int, const std::vector< Trackster > &) const
void calculateDistanceToHigher(const TILES &, const unsigned int layerId, const std::vector< std::pair< int, int >> &)
EventAuxiliary const & eventAuxiliary() const override
double eta() const
pseudorapidity of cluster centroid
constexpr std::array< uint8_t, layerIndexSize > layer
static std::string const input
bool getData(T &iHolder) const
void assignPCAtoTracksters(std::vector< Trackster > &, const std::vector< reco::CaloCluster > &, const edm::ValueMap< std::pair< float, float >> &, double, bool energyWeight=true)
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
void makeTracksters(const typename PatternRecognitionAlgoBaseT< TILES >::Inputs &input, std::vector< Trackster > &result, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation) override
std::vector< float > features(const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
void dumpTiles(const TILES &) const
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Abs< T >::type abs(const T &t)
std::atomic< tensorflow::GraphDef * > eidGraphDef
double energy() const
cluster energy
int findAndAssignTracksters(const TILES &, const std::vector< std::pair< int, int >> &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void calculateLocalDensity(const TILES &, const unsigned int layerId, const std::vector< std::pair< int, int >> &)
PatternRecognitionbyCLUE3D(const edm::ParameterSet &conf, const CacheBase *cache, edm::ConsumesCollector)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
void dumpClusters(const std::vector< std::pair< int, int >> &layerIdx2layerandSoa, const int) const
tensorflow::Session * eidSession_
char data[epos_bytes_allocation]
double phi() const
azimuthal angle of cluster centroid
void reset(double vett[256])
int sum_x
More diagnostics.
EventNumber_t event() const
*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
void energyRegressionAndID(const std::vector< reco::CaloCluster > &layerClusters, std::vector< Trackster > &result)