LDMX Software
CLUE.cxx
1
2#include "Ecal/CLUE.h"
3
4#include <cmath>
5
6namespace ecal {
7
14template <typename T>
15T CLUE::dist(T x1, T y1, T x2, T y2) {
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);
20}
21
22// 3D version, overloaded
23template <typename T>
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);
30}
31
32/* Old code, idea was to do electron reclustering based on first layer
33 centroids' distances to each other I.e. if electrons are close together =>
34 likely merged => recluster Did not quite work and I don't remember the idea
35 anymore but leaving the code here for inspo */
36
37void CLUE::electronSeparation(std::vector<ldmx::EcalHit> hits) {
38 std::vector<double> layer_thickness = {2., 3.5, 5.3, 5.3, 5.3, 5.3,
39 5.3, 5.3, 5.3, 5.3, 5.3, 10.5,
40 10.5, 10.5, 10.5, 10.5};
41 double air = 10.;
42 // sort hits in z
43 std::sort(hits.begin(), hits.end(),
44 [](const ldmx::EcalHit& a, const ldmx::EcalHit& b) {
45 return a.getZPos() < b.getZPos();
46 });
47
48 std::vector<ldmx::EcalHit> first_layers;
49 std::vector<IntermediateCluster> first_layer_clusters;
50 int layer_tag = 0;
51 double layer_z = hits[0].getZPos();
52 for (const auto& hit : hits) {
53 if (hit.getZPos() > layer_z + layer_thickness[layer_tag] + air) {
54 layer_tag++;
55 // if (layerTag > limit) break;
56 break;
57 }
58 first_layers.push_back(hit);
59 first_layer_clusters.push_back(IntermediateCluster(hit, layer_tag));
60 }
61 bool merge = false;
62 do {
63 merge = false;
64 for (int i = 0; i < first_layer_clusters.size(); i++) {
65 if (first_layer_clusters[i].empty()) continue;
66 // if (firstLayerClusters[i].centroid().E() >= seedThreshold_) {
67 for (int j = i + 1; j < first_layer_clusters.size(); j++) {
68 if (first_layer_clusters[j].empty()) continue;
69 if (dist(first_layer_clusters[i].centroid().Px(),
70 first_layer_clusters[i].centroid().Py(),
71 first_layer_clusters[j].centroid().Px(),
72 first_layer_clusters[j].centroid().Py()) < 8.) {
73 first_layer_clusters[i].add(first_layer_clusters[j]);
74 first_layer_clusters[j].clear();
75 merge = true;
76 }
77 }
78 // } else break;
79 }
80 } while (merge);
81 ldmx_log(trace) << "--- ELECTRON SEPARATION ---";
82 for (int i = 0; i < first_layer_clusters.size(); i++) {
83 if (first_layer_clusters[i].empty()) continue;
84 ldmx_log(trace) << " Cluster " << i
85 << " x: " << first_layer_clusters[i].centroid().Px()
86 << " y: " << first_layer_clusters[i].centroid().Py();
87 for (int j = i + 1; j < first_layer_clusters.size(); j++) {
88 if (first_layer_clusters[j].empty()) continue;
89 auto d = dist(first_layer_clusters[i].centroid().Px(),
90 first_layer_clusters[i].centroid().Py(),
91 first_layer_clusters[j].centroid().Px(),
92 first_layer_clusters[j].centroid().Py());
93 ldmx_log(trace) << "Dist to cluster " << j << ": " << d;
94 }
95 }
96}
97
98// Function to create layers from hits based on their layer number
99std::vector<std::vector<const ldmx::EcalHit*>> CLUE::createLayers(
100 const std::vector<const ldmx::EcalHit*>& hits) {
101 ldmx_log(trace) << "--- LAYER CREATION ---";
102 ldmx_log(trace) << "Number of layers: " << nbr_of_layers_;
103
104 // vector of layers, each layer having a vector of hits
105 // initialize with nbr_of_layers_ empty vectors
106 std::vector<std::vector<const ldmx::EcalHit*>> hits_per_layer(nbr_of_layers_);
107
108 // Track highest energy per layer for proper rho_c calculation
109 std::vector<double> layer_max_energies(nbr_of_layers_, 0.0);
110
111 // Clear any existing layer_rho_c_ values
112 layer_rho_c_.clear();
113
114 // This is rather ad-hoc, but this way it's still tunable via rhoc_ parameter
115 double rhoc_factor = rhoc_ / 250.;
116
117 for (const auto& hit : hits) {
118 ldmx::EcalID ecal_id(hit->getID());
119 // layer number from EcalID, starting at 0
120 int layer = ecal_id.layer();
121
122 if (layer >= nbr_of_layers_) {
123 ldmx_log(trace) << "Skipping hit in layer " << layer
124 << " (beyond nbr_of_layers_ = " << nbr_of_layers_ << ")";
125 continue;
126 }
127
128 // Track the highest energy in this specific layer
129 if (hit->getEnergy() > layer_max_energies[layer]) {
130 layer_max_energies[layer] = hit->getEnergy();
131 }
132
133 ldmx_log(trace) << " Adding hit with energy " << hit->getEnergy()
134 << " to layer " << layer;
135 hits_per_layer[layer].push_back(hit);
136 } // end of loop over hits
137
138 // Calculate rhoc for each layer based on that layer's maximum energy
139 for (int i = 0; i < nbr_of_layers_; i++) {
140 double rho_c = layer_max_energies[i] / rhoc_factor;
141 layer_rho_c_.push_back(rho_c);
142 ldmx_log(trace) << "Layer " << i << " has " << hits_per_layer[i].size()
143 << " hits, max energy " << layer_max_energies[i]
144 << ", rho_c " << rho_c;
145 }
146
147 ldmx_log(trace) << "Created " << hits_per_layer.size() << " layers";
148 return hits_per_layer;
149} // end of createLayers
150
151float CLUE::roundToDecimal(float x_, int num_decimal_precision_digits) {
152 float power_of_10 = std::pow(10, num_decimal_precision_digits);
153 return std::round(x_ * power_of_10) / power_of_10;
154}
155
156std::vector<std::shared_ptr<CLUE::Density>> CLUE::setup(
157 const std::vector<const ldmx::EcalHit*>& hits) {
158 std::vector<std::shared_ptr<Density>> densities;
159 std::map<std::pair<float, float>, std::shared_ptr<Density>> density_map;
160 event_centroid_ = IntermediateCluster();
161 ldmx_log(trace) << "--- SETUP ---";
162 ldmx_log(trace) << "Building densities";
163 for (const auto& hit : hits) {
164 // collapse z dimension
165 float x = roundToDecimal(hit->getXPos(), 4);
166 float y = roundToDecimal(hit->getYPos(), 4);
167 float z = roundToDecimal(hit->getZPos(), 4);
168 ldmx_log(trace) << " New hit { x: " << x << " y: " << y << "}"
169 << " (and z: " << z << ")";
170 std::pair<float, float> coords;
171 if (dc_ != 0 && nbr_of_layers_ > 1) {
172 // if more than one layer, divide hit into densities with side dc
173 double i = std::ceil(std::abs(x) / dc_);
174 double j = std::ceil(std::abs(y) / dc_);
175 if (x < 0) {
176 i = -i;
177 x = (i + 0.5) * dc_;
178 } else {
179 x = (i - 0.5) * dc_;
180 }
181 if (y < 0) {
182 j = -j;
183 y = (j + 0.5) * dc_;
184 } else {
185 y = (j - 0.5) * dc_;
186 }
187 coords = {i, j};
188 ldmx_log(trace) << " Index " << i << ", " << j << "; x: " << x
189 << " y: " << y;
190 } else {
191 // if just one layer, have all densities with the same x,y be in same
192 // density
193 coords = {x, y};
194 }
195
196 if (density_map.find(coords) == density_map.end()) {
197 density_map.emplace(coords, std::make_shared<CLUE::Density>(x, y));
198 ldmx_log(trace) << " * New density created";
199 } else {
200 ldmx_log(trace) << " --> Found density with x: "
201 << density_map[coords]->x_
202 << " y: " << density_map[coords]->y_;
203 }
204 density_map[coords]->hits_.push_back(hit);
205 density_map[coords]->total_energy_ += hit->getEnergy();
206 density_map[coords]->z_ += hit->getZPos();
207
208 event_centroid_.add(hit);
209 } // end of loop over hits
210
211 densities.reserve(density_map.size());
212 for (const auto& entry : density_map) {
213 densities.push_back(std::move(entry.second));
214 }
215 // sort according to energy
216 std::sort(densities.begin(), densities.end(),
217 [](const std::shared_ptr<CLUE::Density>& a,
218 const std::shared_ptr<CLUE::Density>& b) {
219 return a->total_energy_ > b->total_energy_;
220 });
221
222 ldmx_log(trace) << "Decide parents";
223
224 // decide delta_ and follower_of_
225 for (int i = 0; i < densities.size(); i++) {
226 densities[i]->index_ = i;
227 // avg z position
228 densities[i]->z_ = densities[i]->z_ / densities[i]->hits_.size();
229 ldmx_log(trace) << " Index: " << i << "; x: " << densities[i]->x_
230 << "; y: " << densities[i]->y_
231 << "; Energy: " << densities[i]->total_energy_;
232 // loop through all higher energy densities
233 for (int j = 0; j < i; j++) {
234 float distance_2d = dist(densities[i]->x_, densities[i]->y_,
235 densities[j]->x_, densities[j]->y_);
236 // condition energyJ > energyI but this should be baked in as we sorted
237 // according to energy
238 if ((distance_2d < dm_) && (distance_2d < densities[i]->delta_)) {
239 densities[i]->delta_ = distance_2d;
240 densities[i]->follower_of_ = j;
241 ldmx_log(trace) << " New parent, index " << j
242 << "; delta_2d: " << std::setprecision(4)
243 << distance_2d;
244 }
245 }
246 }
247 return densities;
248} // end of setup
249
250// connecting_layers marks if we're currently doing 3D clustering (i.e.
251// connecting seeds between layers) otherwise, layer_index tells us which layer
252// number we're working on
253std::vector<std::vector<const ldmx::EcalHit*>> CLUE::clustering(
254 std::vector<std::shared_ptr<CLUE::Density>>& densities,
255 bool connecting_layers, int layer_index) {
256 ldmx_log(trace) << "--- CLUSTERING ---";
257 ldmx_log(trace) << "Number of densities: " << densities.size()
258 << "; connecting_layers: " << connecting_layers
259 << "; layer_index: " << layer_index;
260 if (!connecting_layers && nbr_of_layers_ > 1) {
261 // if layerwise clustering override rhoc_ and deltac_ with per-layer values
262 rhoc_ = layer_rho_c_[layer_index];
263 ldmx_log(trace) << "Setting rho_c on layer " << layer_index << " to "
264 << rhoc_;
265 if (layer_index * 2 - 1 < radius_.size()) {
266 deltac_ = radius_[layer_index * 2 - 1];
267 ldmx_log(trace) << "Setting delta_c on layer " << layer_index << " to "
268 << deltac_;
269 }
270 } else if (connecting_layers) {
271 // if doing 3D clustering, override deltac_ and rhoc_ with other values
272 // NOT IMPLEMENTED YET
273 // if currently doing 3D clustering
274 // IF 3D clustering
275 // deltao_ = 200.;
276 // deltac_ = 100.;
277 // rhoc_ = 1000.;
278 }
279
280 bool energy_overload = false;
281 double max_energy = 10000.;
282 clustering_loops_ = 0;
283 double delta_c_mod = deltac_;
284 double centroid_radius = 10.;
285
286 // stores seeds of this layer
287 std::vector<std::shared_ptr<Density>>& layer_seeds = seeds_[layer_index];
288
289 // stores hits in cluster
290 std::vector<std::vector<const ldmx::EcalHit*>> clusters;
291 // keeps track of which densities have merged; only used if reclustering
292 std::vector<bool> merged_densities; // index_= cluster id
293 merged_densities.resize(densities.size());
294 // keeps track of cluster energies
295 std::vector<double> cluster_energies;
296 do {
297 // while no cluster has merged
298 if (energy_overload) {
299 // makes delta_c smaller if clusters have merged
300 delta_c_mod = delta_c_mod / 1.1;
301 ldmx_log(trace) << "Energy overload, new delta_cmod: " << delta_c_mod;
302 energy_overload = false;
303 }
304
305 clustering_loops_++;
306 ldmx_log(trace) << "Clustering loop " << clustering_loops_;
307
308 // cluster index
309 int k = 0;
310
311 layer_seeds.clear();
312 layer_seeds.reserve(densities.size());
313 clusters.clear();
314 clusters.reserve(densities.size());
315 cluster_energies.clear();
316 cluster_energies.reserve(densities.size());
317
318 std::stack<int> cluster_stack;
319 // stores followers of densities at corr index_
320 std::vector<std::vector<int>> followers;
321 followers.resize(densities.size());
322
323 // Mark as seed, follower, or outlier
324 for (auto& density : densities) {
325 // funky line to generalize this function for both 2D and 3D case
326 ldmx_log(trace) << " Index: " << density->index_
327 << "; x: " << density->x_ << "; y: " << density->y_
328 << "; Energy: " << density->total_energy_
329 << " Parent ID: " << density->follower_of_
330 << "; Delta: " << density->delta_;
331
332 bool is_seed;
333 if (delta_c_mod != deltac_ && merged_densities[density->cluster_id_] &&
334 dist(density->x_, density->y_, event_centroid_.centroid().x(),
335 event_centroid_.centroid().y()) < centroid_radius) {
336 // if energy has been overloaded and this density belongs to cluster
337 // that was overloaded and this density is close enough to event
338 // centroid use modded delta c
339 is_seed =
340 density->total_energy_ > rhoc_ && density->delta_ > delta_c_mod;
341 } else {
342 is_seed = density->total_energy_ > rhoc_ && density->delta_ > deltac_;
343 if (is_seed) {
344 ldmx_log(trace) << " Distance to event centroid: "
345 << dist(density->x_, density->y_,
346 event_centroid_.centroid().x(),
347 event_centroid_.centroid().y());
348 }
349 }
350 bool is_outlier =
351 (density->total_energy_ < rhoc_) && (density->delta_ > deltao_);
352 density->cluster_id_ = -1;
353 if (is_seed) {
354 ldmx_log(trace) << " This is a Seed";
355 ldmx_log(trace) << " Distance to centroid: "
356 << dist(density->x_, density->y_,
357 event_centroid_.centroid().x(),
358 event_centroid_.centroid().y())
359 << "; with delta " << density->delta_;
360 ldmx_log(trace) << " Setting cluster ID to " << k;
361 density->cluster_id_ = k;
362 k++;
363 // get the index of the seed density
364 cluster_stack.push(density->index_);
365 clusters.push_back(density->hits_);
366 cluster_energies.push_back(density->total_energy_);
367 layer_seeds.push_back(density);
368 } else if (!is_outlier) {
369 ldmx_log(trace) << " This is a Follower";
370 int& parent_index = density->follower_of_;
371 if (parent_index != -1)
372 followers[parent_index].push_back(density->index_);
373 else
374 ldmx_log(error)
375 << " Somehow found a follower with parent index -1: id = "
376 << density->index_;
377 } else {
378 ldmx_log(trace) << " This is an Outlier";
379 }
380 }
381
382 merged_densities.clear();
383 merged_densities.resize(densities.size());
384
385 // Go through all seeds and add followers, then follower's followers, etc.
386 while (cluster_stack.size() > 0) {
387 auto& d = densities[cluster_stack.top()];
388 cluster_stack.pop();
389 auto& cid = d->cluster_id_;
390 // for indices of followers of dp
391 for (const auto follower_index : followers[d->index_]) {
392 auto& follower = densities[follower_index];
393 // Set cluster index of the follower to the cluster index of `d`
394 follower->cluster_id_ = cid;
395 cluster_energies[cid] += follower->total_energy_;
396
397 if (reclustering_ && cluster_energies[cid] > max_energy &&
398 delta_c_mod > 0.5 && clustering_loops_ < 100) {
399 // If reclustering is on and cluster energy is too high,
400 // delta_c_mod is not too low, and we haven't tried for too long
401 merged_densities[cid] = true;
402 if (!energy_overload && clustering_loops_ == 99) {
403 ldmx_log(warn) << "Merging clusters, max cluster loops hit";
404 }
405 energy_overload = true;
406 if (clustering_loops_ != 1) {
407 goto endwhile; // Don't break on the first loop to save the initial
408 // cluster number
409 }
410 }
411
412 clusters[cid].insert(std::end(clusters[cid]),
413 std::begin(follower->hits_),
414 std::end(follower->hits_));
415 // Add follower to the stack so its followers can also get the correct
416 // cluster index
417 cluster_stack.push(follower_index);
418 }
419 }
420 // for first clusteringloop, we want to save number of clusters before
421 // reclustering
422 if (clustering_loops_ == 1 && energy_overload)
423 initial_cluster_nbr_ = clusters.size();
424 endwhile:;
425 } while (energy_overload);
426 // if we have more than one layer and we are not currently doing CLUE3D
427 if (!connecting_layers && nbr_of_layers_ > 1) {
428 // Overwrite seed densities' properties with cluster properties
429 // Might be cleaner to just create new densities for cluster seeds
430 for (auto& seed : layer_seeds) {
431 seed->delta_ = std::numeric_limits<float>::max();
432 seed->hits_ = clusters[seed->cluster_id_];
433 seed->total_energy_ = cluster_energies[seed->cluster_id_];
434 seed->index_ = seed_index_;
435 seed_index_++;
436 }
437 // Sort seeds in layer based on energy
438 std::sort(layer_seeds.begin(), layer_seeds.end(),
439 [](const std::shared_ptr<Density>& a,
440 const std::shared_ptr<Density>& b) {
441 return a->total_energy_ > b->total_energy_;
442 });
443 }
444 return clusters;
445} // end of clustering
446
447// Function to connect seeds between layers to form 3D clusters
448// ONLY used in CLUE3D
449std::vector<std::shared_ptr<CLUE::Density>> CLUE::setupForClue3D() {
450 ldmx_log(trace) << "--- LAYER SETUP ---";
451 std::vector<std::shared_ptr<CLUE::Density>> densities;
452 layer_rho_c_.clear();
453 for (int layer = 0; layer < nbr_of_layers_; layer++) {
454 ldmx_log(trace) << " LAYER " << layer << " with " << seeds_[layer].size()
455 << " seeds";
456 auto& seeds_in_current_layer = seeds_[layer];
457 double highest_energy = 0.;
458 for (const auto& current_seed : seeds_in_current_layer) {
459 // for each seed in layer
460 current_seed->layer_ = layer;
461 if (current_seed->total_energy_ > highest_energy)
462 highest_energy = current_seed->total_energy_;
463 ldmx_log(trace) << " Density with index " << current_seed->index_
464 << ", energy: " << current_seed->total_energy_
465 << " position (x,y)= {" << current_seed->x_ << ","
466 << current_seed->y_ << ")";
467 int depth = 1;
468 // decide delta_ and followerof from seeds in previous and next layer_
469 // do {
470 // depth++;
471 if ((layer - depth >= 0) && (layer - depth < seeds_.size())) {
472 ldmx_log(trace) << " Looking at pre-layer: " << layer - depth;
473 // look at previous layer
474 ldmx_log(trace) << " In previous layer... ";
475 auto& previous_layer = seeds_[layer - depth];
476 ldmx_log(trace) << " Got " << previous_layer.size() << " seeds";
477 for (const auto& prev_seed : previous_layer) {
478 // for each seed in previous layer
479 auto distance_2d_prev = dist(current_seed->x_, current_seed->y_,
480 prev_seed->x_, prev_seed->y_);
481 auto dz_prev = std::abs(current_seed->z_ - prev_seed->z_);
482 ldmx_log(trace) << " DeltaXY to index " << prev_seed->index_
483 << ": " << std::setprecision(4) << distance_2d_prev;
484 ldmx_log(trace) << " DeltaZ to index " << prev_seed->index_
485 << ": " << std::setprecision(4) << dz_prev;
486 if (prev_seed->total_energy_ > current_seed->total_energy_ &&
487 distance_2d_prev < current_seed->delta_ &&
488 dz_prev < current_seed->z_delta_) {
489 ldmx_log(trace) << " New parent: index " << prev_seed->index_
490 << " on layer " << layer - depth << "; energy "
491 << prev_seed->total_energy_;
492 ldmx_log(trace) << " New 2D distance to prev-layer: "
493 << std::setprecision(4) << distance_2d_prev;
494 ldmx_log(trace)
495 << " New delta Z: " << std::setprecision(4) << dz_prev;
496 current_seed->delta_ = distance_2d_prev;
497 current_seed->z_delta_ = dz_prev;
498 current_seed->follower_of_ = prev_seed->index_;
499 } else if (prev_seed->total_energy_ < current_seed->total_energy_) {
500 ldmx_log(trace) << " Breaking on index " << prev_seed->index_
501 << " with energy " << prev_seed->total_energy_;
502 break;
503 }
504 }
505 }
506
507 if (layer + depth < nbr_of_layers_ && layer + depth < seeds_.size()) {
508 ldmx_log(trace) << " Looking at post-layer: " << layer + depth;
509 auto& next_layer = seeds_[layer + depth];
510 ldmx_log(trace) << " Got " << next_layer.size() << " seeds";
511 for (const auto& next_seed : next_layer) {
512 auto distance_2d_next = dist(current_seed->x_, current_seed->y_,
513 next_seed->x_, next_seed->y_);
514 auto dz_next = std::abs(current_seed->z_ - next_seed->z_);
515 ldmx_log(trace) << " DeltaXY to index " << next_seed->index_
516 << ": " << std::setprecision(4) << distance_2d_next;
517 ldmx_log(trace) << " DeltaZ to index_" << next_seed->index_
518 << ": " << std::setprecision(4) << dz_next;
519 if (next_seed->total_energy_ > current_seed->total_energy_ &&
520 distance_2d_next < current_seed->delta_ &&
521 dz_next < current_seed->z_delta_) {
522 ldmx_log(trace) << " New parent: index_" << next_seed->index_
523 << " on layer " << layer + depth << "; energy "
524 << next_seed->total_energy_;
525 ldmx_log(trace) << " New 2D distance to next-layer: "
526 << std::setprecision(4) << distance_2d_next;
527 ldmx_log(trace)
528 << " New delta_Z: " << std::setprecision(4) << dz_next;
529 current_seed->delta_ = distance_2d_next;
530 current_seed->z_delta_ = dz_next;
531 current_seed->follower_of_ = next_seed->index_;
532 } else if (next_seed->total_energy_ < current_seed->total_energy_) {
533 ldmx_log(trace) << " Breaking on index_" << next_seed->index_
534 << " with energy " << next_seed->total_energy_;
535 break;
536 }
537 } // end loop on seeds in next layer
538 } // end of looking at next layer
539 // } while (currentLayer[i]->layerFollowerOf == -1 && (layer - depth >=
540 // 0 || layer + depth < nbr_of_layers_));
541 ldmx_log(trace) << " Done setting parents.";
542 densities.push_back(current_seed);
543 }
544 // TODO: This 2 needs to be configurable
545 layer_rho_c_.push_back(highest_energy / 2);
546 } // end loop over layers
547 return densities;
548} // end of setupForClue3D
549
550void CLUE::convertToIntermediateClusters(
551 std::vector<std::vector<const ldmx::EcalHit*>>& clusters) {
552 // Convert to workingecalclusters to ensure compatibility with
553 // EcalClusterProducer
554 for (const auto& cluster : clusters) {
555 IntermediateCluster intermediate_cluster{},
556 intermediate_cluster_first_layer{};
557
558 for (const auto& hit : cluster) {
559 intermediate_cluster.add(hit);
560 // if hit is in first layer, add to first layer cluster
561 ldmx::EcalID ecal_id(hit->getID());
562 auto layer = ecal_id.layer();
563 intermediate_cluster.setLayer(layer);
564 if (layer == 0) {
565 intermediate_cluster_first_layer.add(hit);
566 intermediate_cluster_first_layer.setLayer(layer);
567 }
568 }
569 final_clusters_.push_back(intermediate_cluster);
570 first_layer_centroids_.push_back(intermediate_cluster_first_layer);
571 auto cent_x = intermediate_cluster.centroid().x();
572 auto cent_y = intermediate_cluster.centroid().y();
573 auto event_cent_x = event_centroid_.centroid().x();
574 auto event_cent_y = event_centroid_.centroid().y();
575 const auto& distance = dist(cent_x, cent_y, event_cent_x, event_cent_y);
576 centroid_distances_.push_back(distance);
577 }
578}
579
580void CLUE::cluster(const std::vector<ldmx::EcalHit>& unsorted_hits, double dc,
581 double rc, double delta_c, double delta_o, int nbr_of_layers,
582 bool reclustering) {
583 ldmx_log(info) << "Starting CLUE clustering with parameters:" << "dc " << dc
584 << ", rc " << rc << ", delta_c " << delta_c << ", delta_o "
585 << delta_o << ", nbr_of_layers " << nbr_of_layers
586 << ", reclustering " << reclustering;
587 // cutoff distance for local density
588 dc_ = dc;
589 // min density to promote as seed/max density to demote as outlier
590 rhoc_ = rc;
591 // min separation distance for seeds
592 deltac_ = delta_c;
593 // min separation distance for outliers
594 deltao_ = delta_o;
595 // max distance to parent for both seeds and followers
596 dm_ = std::max(delta_c, delta_o);
597
598 // Recluster merged clusters or not
599 reclustering_ = reclustering;
600 nbr_of_layers_ = nbr_of_layers;
601
602 if (nbr_of_layers_ < 1) {
603 // anything below 1 => include all layers
604 nbr_of_layers_ = max_layers_;
605 } else if (nbr_of_layers_ > max_layers_) {
606 ldmx_log(warn) << "nbr_of_layers_ " << nbr_of_layers_
607 << " exceeds max layers " << max_layers_
608 << ", setting to max layers";
609 nbr_of_layers_ = max_layers_;
610 }
611
612 // first copy *addresses* so we are only ever passing around pointers
613 std::vector<const ldmx::EcalHit*> hits;
614 hits.reserve(unsorted_hits.size());
615 ldmx_log(debug) << "Clustering " << unsorted_hits.size() << " hits";
616 for (const auto& unsorted_hit : unsorted_hits) {
617 hits.push_back(&unsorted_hit);
618 }
619 // sort hits by Z position
620 ldmx_log(debug) << "Sorting hits by Z position";
621 std::sort(hits.begin(), hits.end(),
622 [](const ldmx::EcalHit* a, const ldmx::EcalHit* b) {
623 return a->getZPos() < b->getZPos();
624 });
625
626 seeds_.resize(nbr_of_layers_);
627
628 if (nbr_of_layers_ > 1) {
629 ldmx_log(debug) << "Creating layers";
630 // returns a vector of layers, with each layer having a vector of hits
631 const auto layers = createLayers(hits);
632 ldmx_log(debug) << "Doing layer-wise clustering on " << layers.size()
633 << " layers";
634 for (int i = 0; i < layers.size(); i++) {
635 ldmx_log(trace) << "--- LAYER " << i + 1 << " ---";
636 auto densities = setup(layers[i]);
637 auto clusters = clustering(densities, false, i);
638 convertToIntermediateClusters(clusters);
639 // clustering without 3D
640 }
641 // Below for CLUE3D, comment for just layer clustering
642 // This does not work properly yet
643 // auto densities = setupForClue3D();
644 // auto clusters = clustering(densities, true);
645 // convertToIntermediateClusters(clusters);
646 } else {
647 ldmx_log(debug) << "Only one layer, doing 2D clustering";
648 auto densities = setup(hits);
649 auto clusters = clustering(densities, false);
650 convertToIntermediateClusters(clusters);
651 }
652}
653
654} // namespace ecal
A version of CLUE (CMS) for clustering in ECal.
T dist(T x1, T y1, T x2, T y2)
Euclidean distance between two points.
Definition CLUE.cxx:15
Stores reconstructed hit information from the ECAL.
Definition EcalHit.h:19
Extension of DetectorID providing access to ECal layers and cell numbers in a hex grid.
Definition EcalID.h:20