16 auto delta_x = x1 - x2;
17 auto delta_y = y1 - y2;
18 auto r_square = delta_x * delta_x + delta_y * delta_y;
19 return std::sqrt(r_square);
24T
CLUE::dist(T x1, T y1, T z1, T x2, T y2, T z2) {
25 auto delta_x = x1 - x2;
26 auto delta_y = y1 - y2;
27 auto delta_z = z1 - z2;
28 auto r_square = delta_x * delta_x + delta_y * delta_y + delta_z * delta_z;
29 return std::sqrt(r_square);
97std::vector<std::vector<const ldmx::EcalHit*>> CLUE::createLayers(
98 const std::vector<const ldmx::EcalHit*>& hits) {
99 ldmx_log(trace) <<
"--- LAYER CREATION ---";
100 ldmx_log(trace) <<
"Number of layers: " << nbr_of_layers_;
104 std::vector<std::vector<const ldmx::EcalHit*>> hits_per_layer(nbr_of_layers_);
107 std::vector<double> layer_max_energies(nbr_of_layers_, 0.0);
110 layer_rho_c_.clear();
113 double rhoc_factor = rhoc_ / 250.;
115 for (
const auto& hit : hits) {
118 int layer = ecal_id.layer();
120 if (layer >= nbr_of_layers_) {
121 ldmx_log(trace) <<
"Skipping hit in layer " << layer
122 <<
" (beyond nbr_of_layers_ = " << nbr_of_layers_ <<
")";
127 if (hit->getEnergy() > layer_max_energies[layer]) {
128 layer_max_energies[layer] = hit->getEnergy();
131 ldmx_log(trace) <<
" Adding hit with energy " << hit->getEnergy()
132 <<
" to layer " << layer;
133 hits_per_layer[layer].push_back(hit);
137 for (
int i = 0; i < nbr_of_layers_; i++) {
138 double rho_c = layer_max_energies[i] / rhoc_factor;
139 layer_rho_c_.push_back(rho_c);
140 ldmx_log(trace) <<
"Layer " << i <<
" has " << hits_per_layer[i].size()
141 <<
" hits, max energy " << layer_max_energies[i]
142 <<
", rho_c " << rho_c;
145 ldmx_log(trace) <<
"Created " << hits_per_layer.size() <<
" layers";
146 return hits_per_layer;
149float CLUE::roundToDecimal(
float x_,
int num_decimal_precision_digits) {
150 float power_of_10 = std::pow(10, num_decimal_precision_digits);
151 return std::round(x_ * power_of_10) / power_of_10;
154std::vector<std::shared_ptr<CLUE::Density>> CLUE::setup(
155 const std::vector<const ldmx::EcalHit*>& hits) {
156 std::vector<std::shared_ptr<Density>> densities;
157 std::map<std::pair<float, float>, std::shared_ptr<Density>> density_map;
158 event_centroid_ = IntermediateCluster();
159 ldmx_log(trace) <<
"--- SETUP ---";
160 ldmx_log(trace) <<
"Building densities";
161 for (
const auto& hit : hits) {
163 float x = roundToDecimal(hit->getXPos(), 4);
164 float y = roundToDecimal(hit->getYPos(), 4);
165 ldmx_log(trace) <<
" New hit { x: " << x <<
" y: " << y <<
"}";
166 std::pair<float, float> coords;
167 if (dc_ != 0 && nbr_of_layers_ > 1) {
169 double i = std::ceil(std::abs(x) / dc_);
170 double j = std::ceil(std::abs(y) / dc_);
184 ldmx_log(trace) <<
" Index " << i <<
", " << j <<
"; x: " << x
191 if (density_map.find(coords) == density_map.end()) {
192 density_map.emplace(coords, std::make_shared<CLUE::Density>(x, y));
193 ldmx_log(trace) <<
" New density created";
195 ldmx_log(trace) <<
" Found density with x: " << density_map[coords]->x_
196 <<
" y: " << density_map[coords]->y_;
198 density_map[coords]->hits_.push_back(hit);
199 density_map[coords]->total_energy_ += hit->getEnergy();
200 density_map[coords]->z_ += hit->getZPos();
202 event_centroid_.add(hit);
205 densities.reserve(density_map.size());
206 for (
const auto& entry : density_map) {
207 densities.push_back(std::move(entry.second));
210 std::sort(densities.begin(), densities.end(),
211 [](
const std::shared_ptr<CLUE::Density>& a,
212 const std::shared_ptr<CLUE::Density>& b) {
213 return a->total_energy_ > b->total_energy_;
216 ldmx_log(trace) <<
"Decide parents";
219 for (
int i = 0; i < densities.size(); i++) {
220 densities[i]->index_ = i;
222 densities[i]->z_ = densities[i]->z_ / densities[i]->hits_.size();
223 ldmx_log(trace) <<
" Index: " << i <<
"; x: " << densities[i]->x_
224 <<
"; y: " << densities[i]->y_
225 <<
"; Energy: " << densities[i]->total_energy_;
227 for (
int j = 0; j < i; j++) {
228 float distance_2d =
dist(densities[i]->x_, densities[i]->y_,
229 densities[j]->x_, densities[j]->y_);
232 if ((distance_2d < dm_) && (distance_2d < densities[i]->delta_)) {
233 densities[i]->delta_ = distance_2d;
234 densities[i]->follower_of_ = j;
235 ldmx_log(trace) <<
" New parent, index " << j
236 <<
"; delta_2d: " << std::setprecision(4)
247std::vector<std::vector<const ldmx::EcalHit*>> CLUE::clustering(
248 std::vector<std::shared_ptr<CLUE::Density>>& densities,
249 bool connecting_layers,
int layer_index) {
250 ldmx_log(trace) <<
"--- CLUSTERING ---";
251 ldmx_log(trace) <<
"Number of densities: " << densities.size()
252 <<
"; connecting_layers: " << connecting_layers
253 <<
"; layer_index: " << layer_index;
254 if (!connecting_layers && nbr_of_layers_ > 1) {
256 rhoc_ = layer_rho_c_[layer_index];
257 ldmx_log(trace) <<
"Setting rho_c on layer " << layer_index <<
" to "
259 if (layer_index * 2 - 1 < radius_.size()) {
260 deltac_ = radius_[layer_index * 2 - 1];
261 ldmx_log(trace) <<
"Setting delta_c on layer " << layer_index <<
" to "
264 }
else if (connecting_layers) {
274 bool energy_overload =
false;
275 double max_energy = 10000.;
276 clustering_loops_ = 0;
277 double delta_c_mod = deltac_;
278 double centroid_radius = 10.;
281 std::vector<std::shared_ptr<Density>>& layer_seeds = seeds_[layer_index];
284 std::vector<std::vector<const ldmx::EcalHit*>> clusters;
286 std::vector<bool> merged_densities;
287 merged_densities.resize(densities.size());
289 std::vector<double> cluster_energies;
292 if (energy_overload) {
294 delta_c_mod = delta_c_mod / 1.1;
295 ldmx_log(trace) <<
"Energy overload, new delta_cmod: " << delta_c_mod;
296 energy_overload =
false;
300 ldmx_log(trace) <<
"Clustering loop " << clustering_loops_;
306 layer_seeds.reserve(densities.size());
308 clusters.reserve(densities.size());
309 cluster_energies.clear();
310 cluster_energies.reserve(densities.size());
312 std::stack<int> cluster_stack;
314 std::vector<std::vector<int>> followers;
315 followers.resize(densities.size());
318 for (
auto& density : densities) {
320 ldmx_log(trace) <<
" Index: " << density->index_
321 <<
"; x: " << density->x_ <<
"; y: " << density->y_
322 <<
"; Energy: " << density->total_energy_
323 <<
" Parent ID: " << density->follower_of_
324 <<
"; Delta: " << density->delta_;
327 if (delta_c_mod != deltac_ && merged_densities[density->cluster_id_] &&
328 dist(density->x_, density->y_, event_centroid_.centroid().x(),
329 event_centroid_.centroid().y()) < centroid_radius) {
334 density->total_energy_ > rhoc_ && density->delta_ > delta_c_mod;
336 is_seed = density->total_energy_ > rhoc_ && density->delta_ > deltac_;
340 (density->total_energy_ < rhoc_) && (density->delta_ > deltao_);
341 density->cluster_id_ = -1;
343 ldmx_log(trace) <<
" This is a Seed";
344 ldmx_log(trace) <<
" Distance to centroid: "
345 <<
dist(density->x_, density->y_,
346 event_centroid_.centroid().x(),
347 event_centroid_.centroid().y())
348 <<
"; with delta " << density->delta_;
349 ldmx_log(trace) <<
" Setting cluster ID to " << k;
350 density->cluster_id_ = k;
353 cluster_stack.push(density->index_);
354 clusters.push_back(density->hits_);
355 cluster_energies.push_back(density->total_energy_);
356 layer_seeds.push_back(density);
357 }
else if (!is_outlier) {
358 ldmx_log(trace) <<
" This is a Follower";
359 int& parent_index = density->follower_of_;
360 if (parent_index != -1)
361 followers[parent_index].push_back(density->index_);
363 ldmx_log(trace) <<
" This is an Outlier";
367 merged_densities.clear();
368 merged_densities.resize(densities.size());
371 while (cluster_stack.size() > 0) {
372 auto& d = densities[cluster_stack.top()];
374 auto& cid = d->cluster_id_;
376 for (
const auto follower_index : followers[d->index_]) {
377 auto& follower = densities[follower_index];
379 follower->cluster_id_ = cid;
380 cluster_energies[cid] += follower->total_energy_;
382 if (reclustering_ && cluster_energies[cid] > max_energy &&
383 delta_c_mod > 0.5 && clustering_loops_ < 100) {
386 merged_densities[cid] =
true;
387 if (!energy_overload && clustering_loops_ == 99) {
388 ldmx_log(warn) <<
"Merging clusters, max cluster loops hit";
390 energy_overload =
true;
391 if (clustering_loops_ != 1) {
397 clusters[cid].insert(std::end(clusters[cid]),
398 std::begin(follower->hits_),
399 std::end(follower->hits_));
402 cluster_stack.push(follower_index);
407 if (clustering_loops_ == 1 && energy_overload)
408 initial_cluster_nbr_ = clusters.size();
410 }
while (energy_overload);
412 if (!connecting_layers && nbr_of_layers_ > 1) {
415 for (
auto& seed : layer_seeds) {
416 seed->delta_ = std::numeric_limits<float>::max();
417 seed->hits_ = clusters[seed->cluster_id_];
418 seed->total_energy_ = cluster_energies[seed->cluster_id_];
419 seed->index_ = seed_index_;
423 std::sort(layer_seeds.begin(), layer_seeds.end(),
424 [](
const std::shared_ptr<Density>& a,
425 const std::shared_ptr<Density>& b) {
426 return a->total_energy_ > b->total_energy_;
434std::vector<std::shared_ptr<CLUE::Density>> CLUE::setupForClue3D() {
435 ldmx_log(trace) <<
"--- LAYER SETUP ---";
436 std::vector<std::shared_ptr<CLUE::Density>> densities;
437 layer_rho_c_.clear();
438 for (
int layer = 0; layer < nbr_of_layers_; layer++) {
439 ldmx_log(trace) <<
" LAYER " << layer <<
" with " << seeds_[layer].size()
441 auto& seeds_in_current_layer = seeds_[layer];
442 double highest_energy = 0.;
443 for (
const auto& current_seed : seeds_in_current_layer) {
445 current_seed->layer_ = layer;
446 if (current_seed->total_energy_ > highest_energy)
447 highest_energy = current_seed->total_energy_;
448 ldmx_log(trace) <<
" Density with index " << current_seed->index_
449 <<
", energy: " << current_seed->total_energy_;
454 if ((layer - depth >= 0) && (layer - depth < seeds_.size())) {
455 ldmx_log(trace) <<
" Looking at pre-layer: " << layer - depth;
457 auto& previous_layer = seeds_[layer - depth];
458 for (
const auto& prev_seed : previous_layer) {
460 auto distance_2d_prev =
dist(current_seed->x_, current_seed->y_,
461 prev_seed->x_, prev_seed->y_);
462 auto dz_prev = std::abs(current_seed->z_ - prev_seed->z_);
463 ldmx_log(trace) <<
" DeltaXY to index " << prev_seed->index_
464 <<
": " << std::setprecision(4) << distance_2d_prev;
465 ldmx_log(trace) <<
" DeltaZ to index " << prev_seed->index_
466 <<
": " << std::setprecision(4) << dz_prev;
467 if (prev_seed->total_energy_ > current_seed->total_energy_ &&
468 distance_2d_prev < current_seed->delta_ &&
469 dz_prev < current_seed->z_delta_) {
470 ldmx_log(trace) <<
" New parent: index " << prev_seed->index_
471 <<
" on layer " << layer - depth <<
"; energy "
472 << prev_seed->total_energy_;
473 ldmx_log(trace) <<
" New 2D distance to prev-layer: "
474 << std::setprecision(4) << distance_2d_prev;
476 <<
" New delta Z: " << std::setprecision(4) << dz_prev;
477 current_seed->delta_ = distance_2d_prev;
478 current_seed->z_delta_ = dz_prev;
479 current_seed->follower_of_ = prev_seed->index_;
480 }
else if (prev_seed->total_energy_ < current_seed->total_energy_) {
481 ldmx_log(trace) <<
" Breaking on index " << prev_seed->index_
482 <<
" with energy " << prev_seed->total_energy_;
487 if (layer + depth < nbr_of_layers_ && layer + depth < seeds_.size()) {
488 ldmx_log(trace) <<
" Looking at post-layer: " << layer + depth;
489 auto& next_layer = seeds_[layer + depth];
490 for (
const auto& next_seed : next_layer) {
491 auto distance_2d_next =
dist(current_seed->x_, current_seed->y_,
492 next_seed->x_, next_seed->y_);
493 auto dz_next = std::abs(current_seed->z_ - next_seed->z_);
494 ldmx_log(trace) <<
" DeltaXY to index " << next_seed->index_
495 <<
": " << std::setprecision(4) << distance_2d_next;
496 ldmx_log(trace) <<
" DeltaZ to index_" << next_seed->index_
497 <<
": " << std::setprecision(4) << dz_next;
498 if (next_seed->total_energy_ > current_seed->total_energy_ &&
499 distance_2d_next < current_seed->delta_ &&
500 dz_next < current_seed->z_delta_) {
501 ldmx_log(trace) <<
" New parent: index_" << next_seed->index_
502 <<
" on layer " << layer + depth <<
"; energy "
503 << next_seed->total_energy_;
504 ldmx_log(trace) <<
" New 2D distance to next-layer: "
505 << std::setprecision(4) << distance_2d_next;
507 <<
" New delta_Z: " << std::setprecision(4) << dz_next;
508 current_seed->delta_ = distance_2d_next;
509 current_seed->z_delta_ = dz_next;
510 current_seed->follower_of_ = next_seed->index_;
511 }
else if (next_seed->total_energy_ < current_seed->total_energy_) {
512 ldmx_log(trace) <<
" Breaking on index_" << next_seed->index_
513 <<
" with energy " << next_seed->total_energy_;
520 densities.push_back(current_seed);
523 layer_rho_c_.push_back(highest_energy / 2);
528void CLUE::convertToIntermediateClusters(
529 std::vector<std::vector<const ldmx::EcalHit*>>& clusters) {
532 for (
const auto& cluster : clusters) {
533 IntermediateCluster intermediate_cluster{},
534 intermediate_cluster_first_layer{};
536 for (
const auto& hit : cluster) {
537 intermediate_cluster.add(hit);
540 auto layer = ecal_id.layer();
541 intermediate_cluster.setLayer(layer);
543 intermediate_cluster_first_layer.add(hit);
544 intermediate_cluster_first_layer.setLayer(layer);
547 final_clusters_.push_back(intermediate_cluster);
548 first_layer_centroids_.push_back(intermediate_cluster_first_layer);
549 auto cent_x = intermediate_cluster.centroid().x();
550 auto cent_y = intermediate_cluster.centroid().y();
551 auto event_cent_x = event_centroid_.centroid().x();
552 auto event_cent_y = event_centroid_.centroid().y();
553 const auto& distance =
dist(cent_x, cent_y, event_cent_x, event_cent_y);
554 centroid_distances_.push_back(distance);
558void CLUE::cluster(
const std::vector<ldmx::EcalHit>& unsorted_hits,
double dc,
559 double rc,
double delta_c,
double delta_o,
int nbr_of_layers,
561 ldmx_log(info) <<
"Starting CLUE clustering with parameters:" <<
"dc " << dc
562 <<
", rc " << rc <<
", delta_c " << delta_c <<
", delta_o "
563 << delta_o <<
", nbr_of_layers " << nbr_of_layers
564 <<
", reclustering " << reclustering;
574 dm_ = std::max(delta_c, delta_o);
577 reclustering_ = reclustering;
578 nbr_of_layers_ = nbr_of_layers;
580 if (nbr_of_layers_ < 1) {
582 nbr_of_layers_ = max_layers_;
583 }
else if (nbr_of_layers_ > max_layers_) {
584 ldmx_log(warn) <<
"nbr_of_layers_ " << nbr_of_layers_
585 <<
" exceeds max layers " << max_layers_
586 <<
", setting to max layers";
587 nbr_of_layers_ = max_layers_;
591 std::vector<const ldmx::EcalHit*> hits;
592 hits.reserve(unsorted_hits.size());
593 ldmx_log(debug) <<
"Clustering " << unsorted_hits.size() <<
" hits";
594 for (
const auto& unsorted_hit : unsorted_hits) {
595 hits.push_back(&unsorted_hit);
598 ldmx_log(debug) <<
"Sorting hits by Z position";
599 std::sort(hits.begin(), hits.end(),
601 return a->getZPos() < b->getZPos();
604 seeds_.resize(nbr_of_layers);
606 if (nbr_of_layers_ > 1) {
607 ldmx_log(debug) <<
"Creating layers";
609 const auto layers = createLayers(hits);
610 ldmx_log(debug) <<
"Doing layer-wise clustering on " << layers.size()
612 for (
int i = 0; i < layers.size(); i++) {
613 ldmx_log(trace) <<
"--- LAYER " << i + 1 <<
" ---";
614 auto densities = setup(layers[i]);
615 auto clusters = clustering(densities,
false, i);
616 convertToIntermediateClusters(clusters);
625 ldmx_log(debug) <<
"Only one layer, doing 2D clustering";
626 auto densities = setup(hits);
627 auto clusters = clustering(densities,
false);
628 convertToIntermediateClusters(clusters);
A version of CLUE (CMS) for clustering in ECal.
T dist(T x1, T y1, T x2, T y2)
Euclidean distance between two points.
Stores reconstructed hit information from the ECAL.
Extension of DetectorID providing access to ECal layers and cell numbers in a hex grid.