LDMX Software
ecal::CLUE Class Reference

Classes

struct  Density
 

Public Member Functions

template<typename T >
dist (T x1, T y1, T x2, T y2)
 Euclidean distance between two points.
 
template<typename T >
dist (T x1, T y1, T z1, T x2, T y2, T z2)
 
std::vector< std::vector< const ldmx::EcalHit * > > createLayers (const std::vector< const ldmx::EcalHit * > &hits)
 
float roundToDecimal (float x, int num_decimal_precision_digits)
 
std::vector< std::shared_ptr< Density > > setup (const std::vector< const ldmx::EcalHit * > &hits)
 
void electronSeparation (std::vector< ldmx::EcalHit > hits)
 
std::vector< std::vector< const ldmx::EcalHit * > > clustering (std::vector< std::shared_ptr< Density > > &densities, bool connectingLayers, int layerTag=0)
 
std::vector< std::shared_ptr< Density > > setupForClue3D ()
 
void convertToIntermediateClusters (std::vector< std::vector< const ldmx::EcalHit * > > &clusters)
 
void cluster (const std::vector< ldmx::EcalHit > &hits, double dc, double rc, double deltac, double deltao, int nbrOfLayers, bool reclustering)
 
std::vector< double > getCentroidDistances () const
 
int getNLoops () const
 
int getInitialClusterNbr () const
 
std::vector< IntermediateClustergetClusters () const
 
std::vector< IntermediateClustergetFirstLayerCentroids () const
 

Private Member Functions

 enableLogging ("CLUE")
 

Private Attributes

int clustering_loops_
 
bool reclustering_
 
double dc_
 
double rhoc_
 
double deltac_
 
double deltao_
 
double dm_
 
int max_layers_ {32}
 
int nbr_of_layers_
 
std::vector< double > layer_rho_c_
 
std::vector< double > layer_delta_c_
 
std::vector< double > radius_
 
std::vector< double > centroid_distances_
 
IntermediateCluster event_centroid_
 
std::vector< IntermediateClusterfirst_layer_centroids_
 
int seed_index_ {0}
 
std::vector< std::vector< std::shared_ptr< Density > > > seeds_
 
int initial_cluster_nbr_ {-1}
 
std::vector< IntermediateClusterfinal_clusters_
 
std::vector< std::pair< double, double > > layer_centroid_separations_
 

Detailed Description

Definition at line 28 of file CLUE.h.

Member Function Documentation

◆ cluster()

void ecal::CLUE::cluster ( const std::vector< ldmx::EcalHit > & hits,
double dc,
double rc,
double deltac,
double deltao,
int nbrOfLayers,
bool reclustering )

Definition at line 582 of file CLUE.cxx.

584 {
585 ldmx_log(info) << "Starting CLUE clustering with parameters:" << "dc " << dc
586 << ", rc " << rc << ", delta_c " << delta_c << ", delta_o "
587 << delta_o << ", nbr_of_layers " << nbr_of_layers
588 << ", reclustering " << reclustering;
589 // cutoff distance for local density
590 dc_ = dc;
591 // min density to promote as seed/max density to demote as outlier
592 rhoc_ = rc;
593 // min separation distance for seeds
594 deltac_ = delta_c;
595 // min separation distance for outliers
596 deltao_ = delta_o;
597 // max distance to parent for both seeds and followers
598 dm_ = std::max(delta_c, delta_o);
599
600 // Recluster merged clusters or not
601 reclustering_ = reclustering;
602 nbr_of_layers_ = nbr_of_layers;
603
604 if (nbr_of_layers_ < 1) {
605 // anything below 1 => include all layers
606 nbr_of_layers_ = max_layers_;
607 } else if (nbr_of_layers_ > max_layers_) {
608 ldmx_log(warn) << "nbr_of_layers_ " << nbr_of_layers_
609 << " exceeds max layers " << max_layers_
610 << ", setting to max layers";
611 nbr_of_layers_ = max_layers_;
612 }
613
614 // first copy *addresses* so we are only ever passing around pointers
615 std::vector<const ldmx::EcalHit*> hits;
616 hits.reserve(unsorted_hits.size());
617 ldmx_log(debug) << "Clustering " << unsorted_hits.size() << " hits";
618 for (const auto& unsorted_hit : unsorted_hits) {
619 hits.push_back(&unsorted_hit);
620 }
621 // sort hits by Z position
622 ldmx_log(debug) << "Sorting hits by Z position";
623 std::sort(hits.begin(), hits.end(),
624 [](const ldmx::EcalHit* a, const ldmx::EcalHit* b) {
625 return a->getZPos() < b->getZPos();
626 });
627
628 seeds_.resize(nbr_of_layers_);
629
630 if (nbr_of_layers_ > 1) {
631 ldmx_log(debug) << "Creating layers";
632 // returns a vector of layers, with each layer having a vector of hits
633 const auto layers = createLayers(hits);
634 ldmx_log(debug) << "Doing layer-wise clustering on " << layers.size()
635 << " layers";
636 for (int i = 0; i < layers.size(); i++) {
637 ldmx_log(trace) << "--- LAYER " << i + 1 << " ---";
638 auto densities = setup(layers[i]);
639 auto clusters = clustering(densities, false, i);
640 convertToIntermediateClusters(clusters);
641 // clustering without 3D
642 }
643 // Below for CLUE3D, comment for just layer clustering
644 // This does not work properly yet
645 // auto densities = setupForClue3D();
646 // auto clusters = clustering(densities, true);
647 // convertToIntermediateClusters(clusters);
648 } else {
649 ldmx_log(debug) << "Only one layer, doing 2D clustering";
650 auto densities = setup(hits);
651 auto clusters = clustering(densities, false);
652 convertToIntermediateClusters(clusters);
653 }
654}
Stores reconstructed hit information from the ECAL.
Definition EcalHit.h:19

◆ clustering()

std::vector< std::vector< const ldmx::EcalHit * > > ecal::CLUE::clustering ( std::vector< std::shared_ptr< Density > > & densities,
bool connectingLayers,
int layerTag = 0 )

Definition at line 255 of file CLUE.cxx.

257 {
258 ldmx_log(trace) << "--- CLUSTERING ---";
259 ldmx_log(trace) << "Number of densities: " << densities.size()
260 << "; connecting_layers: " << connecting_layers
261 << "; layer_index: " << layer_index;
262 if (!connecting_layers && nbr_of_layers_ > 1) {
263 // if layerwise clustering override rhoc_ and deltac_ with per-layer values
264 rhoc_ = layer_rho_c_[layer_index];
265 ldmx_log(trace) << "Setting rho_c on layer " << layer_index << " to "
266 << rhoc_;
267 if (layer_index * 2 - 1 < radius_.size()) {
268 deltac_ = radius_[layer_index * 2 - 1];
269 ldmx_log(trace) << "Setting delta_c on layer " << layer_index << " to "
270 << deltac_;
271 }
272 } else if (connecting_layers) {
273 // if doing 3D clustering, override deltac_ and rhoc_ with other values
274 // NOT IMPLEMENTED YET
275 // if currently doing 3D clustering
276 // IF 3D clustering
277 // deltao_ = 200.;
278 // deltac_ = 100.;
279 // rhoc_ = 1000.;
280 }
281
282 bool energy_overload = false;
283 double max_energy = 10000.;
284 clustering_loops_ = 0;
285 double delta_c_mod = deltac_;
286 double centroid_radius = 10.;
287
288 // stores seeds of this layer
289 std::vector<std::shared_ptr<Density>>& layer_seeds = seeds_[layer_index];
290
291 // stores hits in cluster
292 std::vector<std::vector<const ldmx::EcalHit*>> clusters;
293 // keeps track of which densities have merged; only used if reclustering
294 std::vector<bool> merged_densities; // index_= cluster id
295 merged_densities.resize(densities.size());
296 // keeps track of cluster energies
297 std::vector<double> cluster_energies;
298 do {
299 // while no cluster has merged
300 if (energy_overload) {
301 // makes delta_c smaller if clusters have merged
302 delta_c_mod = delta_c_mod / 1.1;
303 ldmx_log(trace) << "Energy overload, new delta_cmod: " << delta_c_mod;
304 energy_overload = false;
305 }
306
307 clustering_loops_++;
308 ldmx_log(trace) << "Clustering loop " << clustering_loops_;
309
310 // cluster index
311 int k = 0;
312
313 layer_seeds.clear();
314 layer_seeds.reserve(densities.size());
315 clusters.clear();
316 clusters.reserve(densities.size());
317 cluster_energies.clear();
318 cluster_energies.reserve(densities.size());
319
320 std::stack<int> cluster_stack;
321 // stores followers of densities at corr index_
322 std::vector<std::vector<int>> followers;
323 followers.resize(densities.size());
324
325 // Mark as seed, follower, or outlier
326 for (auto& density : densities) {
327 // funky line to generalize this function for both 2D and 3D case
328 ldmx_log(trace) << " Index: " << density->index_
329 << "; x: " << density->x_ << "; y: " << density->y_
330 << "; Energy: " << density->total_energy_
331 << " Parent ID: " << density->follower_of_
332 << "; Delta: " << density->delta_;
333
334 bool is_seed;
335 if (delta_c_mod != deltac_ && merged_densities[density->cluster_id_] &&
336 dist(density->x_, density->y_, event_centroid_.centroidX(),
337 event_centroid_.centroidY()) < centroid_radius) {
338 // if energy has been overloaded and this density belongs to cluster
339 // that was overloaded and this density is close enough to event
340 // centroid use modded delta c
341 is_seed =
342 density->total_energy_ > rhoc_ && density->delta_ > delta_c_mod;
343 } else {
344 is_seed = density->total_energy_ > rhoc_ && density->delta_ > deltac_;
345 if (is_seed) {
346 ldmx_log(trace) << " Distance to event centroid: "
347 << dist(density->x_, density->y_,
348 event_centroid_.centroidX(),
349 event_centroid_.centroidY());
350 }
351 }
352 bool is_outlier =
353 (density->total_energy_ < rhoc_) && (density->delta_ > deltao_);
354 density->cluster_id_ = -1;
355 if (is_seed) {
356 ldmx_log(trace) << " This is a Seed";
357 ldmx_log(trace) << " Distance to centroid: "
358 << dist(density->x_, density->y_,
359 event_centroid_.centroidX(),
360 event_centroid_.centroidY())
361 << "; with delta " << density->delta_;
362 ldmx_log(trace) << " Setting cluster ID to " << k;
363 density->cluster_id_ = k;
364 k++;
365 // get the index of the seed density
366 cluster_stack.push(density->index_);
367 clusters.push_back(density->hits_);
368 cluster_energies.push_back(density->total_energy_);
369 layer_seeds.push_back(density);
370 } else if (!is_outlier) {
371 ldmx_log(trace) << " This is a Follower";
372 int& parent_index = density->follower_of_;
373 if (parent_index != -1)
374 followers[parent_index].push_back(density->index_);
375 else
376 ldmx_log(error)
377 << " Somehow found a follower with parent index -1: id = "
378 << density->index_;
379 } else {
380 ldmx_log(trace) << " This is an Outlier";
381 }
382 }
383
384 merged_densities.clear();
385 merged_densities.resize(densities.size());
386
387 // Go through all seeds and add followers, then follower's followers, etc.
388 while (cluster_stack.size() > 0) {
389 auto& d = densities[cluster_stack.top()];
390 cluster_stack.pop();
391 auto& cid = d->cluster_id_;
392 // for indices of followers of dp
393 for (const auto follower_index : followers[d->index_]) {
394 auto& follower = densities[follower_index];
395 // Set cluster index of the follower to the cluster index of `d`
396 follower->cluster_id_ = cid;
397 cluster_energies[cid] += follower->total_energy_;
398
399 if (reclustering_ && cluster_energies[cid] > max_energy &&
400 delta_c_mod > 0.5 && clustering_loops_ < 100) {
401 // If reclustering is on and cluster energy is too high,
402 // delta_c_mod is not too low, and we haven't tried for too long
403 merged_densities[cid] = true;
404 if (!energy_overload && clustering_loops_ == 99) {
405 ldmx_log(warn) << "Merging clusters, max cluster loops hit";
406 }
407 energy_overload = true;
408 if (clustering_loops_ != 1) {
409 goto endwhile; // Don't break on the first loop to save the initial
410 // cluster number
411 }
412 }
413
414 clusters[cid].insert(std::end(clusters[cid]),
415 std::begin(follower->hits_),
416 std::end(follower->hits_));
417 // Add follower to the stack so its followers can also get the correct
418 // cluster index
419 cluster_stack.push(follower_index);
420 }
421 }
422 // for first clusteringloop, we want to save number of clusters before
423 // reclustering
424 if (clustering_loops_ == 1 && energy_overload)
425 initial_cluster_nbr_ = clusters.size();
426 endwhile:;
427 } while (energy_overload);
428 // if we have more than one layer and we are not currently doing CLUE3D
429 if (!connecting_layers && nbr_of_layers_ > 1) {
430 // Overwrite seed densities' properties with cluster properties
431 // Might be cleaner to just create new densities for cluster seeds
432 for (auto& seed : layer_seeds) {
433 seed->delta_ = std::numeric_limits<float>::max();
434 seed->hits_ = clusters[seed->cluster_id_];
435 seed->total_energy_ = cluster_energies[seed->cluster_id_];
436 seed->index_ = seed_index_;
437 seed_index_++;
438 }
439 // Sort seeds in layer based on energy
440 std::sort(layer_seeds.begin(), layer_seeds.end(),
441 [](const std::shared_ptr<Density>& a,
442 const std::shared_ptr<Density>& b) {
443 return a->total_energy_ > b->total_energy_;
444 });
445 }
446 return clusters;
447} // end of clustering
T dist(T x1, T y1, T x2, T y2)
Euclidean distance between two points.
Definition CLUE.cxx:15
double centroidY() const
Get the centroid Y position (energy-weighted).
double centroidX() const
Get the centroid X position (energy-weighted).

◆ convertToIntermediateClusters()

void ecal::CLUE::convertToIntermediateClusters ( std::vector< std::vector< const ldmx::EcalHit * > > & clusters)

Definition at line 552 of file CLUE.cxx.

553 {
554 // Convert to workingecalclusters to ensure compatibility with
555 // EcalClusterProducer
556 for (const auto& cluster : clusters) {
557 IntermediateCluster intermediate_cluster{},
558 intermediate_cluster_first_layer{};
559
560 for (const auto& hit : cluster) {
561 intermediate_cluster.add(hit);
562 // if hit is in first layer, add to first layer cluster
563 ldmx::EcalID ecal_id(hit->getID());
564 auto layer = ecal_id.layer();
565 intermediate_cluster.setLayer(layer);
566 if (layer == 0) {
567 intermediate_cluster_first_layer.add(hit);
568 intermediate_cluster_first_layer.setLayer(layer);
569 }
570 }
571 final_clusters_.push_back(intermediate_cluster);
572 first_layer_centroids_.push_back(intermediate_cluster_first_layer);
573 auto cent_x = intermediate_cluster.centroidX();
574 auto cent_y = intermediate_cluster.centroidY();
575 auto event_cent_x = event_centroid_.centroidX();
576 auto event_cent_y = event_centroid_.centroidY();
577 const auto& distance = dist(cent_x, cent_y, event_cent_x, event_cent_y);
578 centroid_distances_.push_back(distance);
579 }
580}
recon::WorkingCluster< ldmx::EcalHit > IntermediateCluster
Type alias for WorkingCluster specialized for EcalHit.
Extension of DetectorID providing access to ECal layers and cell numbers in a hex grid.
Definition EcalID.h:20

◆ createLayers()

std::vector< std::vector< const ldmx::EcalHit * > > ecal::CLUE::createLayers ( const std::vector< const ldmx::EcalHit * > & hits)

Definition at line 101 of file CLUE.cxx.

102 {
103 ldmx_log(trace) << "--- LAYER CREATION ---";
104 ldmx_log(trace) << "Number of layers: " << nbr_of_layers_;
105
106 // vector of layers, each layer having a vector of hits
107 // initialize with nbr_of_layers_ empty vectors
108 std::vector<std::vector<const ldmx::EcalHit*>> hits_per_layer(nbr_of_layers_);
109
110 // Track highest energy per layer for proper rho_c calculation
111 std::vector<double> layer_max_energies(nbr_of_layers_, 0.0);
112
113 // Clear any existing layer_rho_c_ values
114 layer_rho_c_.clear();
115
116 // This is rather ad-hoc, but this way it's still tunable via rhoc_ parameter
117 double rhoc_factor = rhoc_ / 250.;
118
119 for (const auto& hit : hits) {
120 ldmx::EcalID ecal_id(hit->getID());
121 // layer number from EcalID, starting at 0
122 int layer = ecal_id.layer();
123
124 if (layer >= nbr_of_layers_) {
125 ldmx_log(trace) << "Skipping hit in layer " << layer
126 << " (beyond nbr_of_layers_ = " << nbr_of_layers_ << ")";
127 continue;
128 }
129
130 // Track the highest energy in this specific layer
131 if (hit->getEnergy() > layer_max_energies[layer]) {
132 layer_max_energies[layer] = hit->getEnergy();
133 }
134
135 ldmx_log(trace) << " Adding hit with energy " << hit->getEnergy()
136 << " to layer " << layer;
137 hits_per_layer[layer].push_back(hit);
138 } // end of loop over hits
139
140 // Calculate rhoc for each layer based on that layer's maximum energy
141 for (int i = 0; i < nbr_of_layers_; i++) {
142 double rho_c = layer_max_energies[i] / rhoc_factor;
143 layer_rho_c_.push_back(rho_c);
144 ldmx_log(trace) << "Layer " << i << " has " << hits_per_layer[i].size()
145 << " hits, max energy " << layer_max_energies[i]
146 << ", rho_c " << rho_c;
147 }
148
149 ldmx_log(trace) << "Created " << hits_per_layer.size() << " layers";
150 return hits_per_layer;
151} // end of createLayers

◆ dist() [1/2]

template<typename T >
T ecal::CLUE::dist ( T x1,
T y1,
T x2,
T y2 )

Euclidean distance between two points.

Euclidean distance between two points Using x*x and std::sqrt since they are more specialized and faster than std::pow.

Template Parameters
Tfloating point type

Definition at line 15 of file CLUE.cxx.

15 {
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}

◆ dist() [2/2]

template<typename T >
T ecal::CLUE::dist ( T x1,
T y1,
T z1,
T x2,
T y2,
T z2 )

Definition at line 24 of file CLUE.cxx.

24 {
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}

◆ electronSeparation()

void ecal::CLUE::electronSeparation ( std::vector< ldmx::EcalHit > hits)

Definition at line 37 of file CLUE.cxx.

37 {
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 IntermediateCluster cluster(hit);
60 cluster.setLayer(layer_tag);
61 first_layer_clusters.push_back(cluster);
62 }
63 bool merge = false;
64 do {
65 merge = false;
66 for (int i = 0; i < first_layer_clusters.size(); i++) {
67 if (first_layer_clusters[i].empty()) continue;
68 // if (firstLayerClusters[i].energy() >= seedThreshold_) {
69 for (int j = i + 1; j < first_layer_clusters.size(); j++) {
70 if (first_layer_clusters[j].empty()) continue;
71 if (dist(first_layer_clusters[i].centroidX(),
72 first_layer_clusters[i].centroidY(),
73 first_layer_clusters[j].centroidX(),
74 first_layer_clusters[j].centroidY()) < 8.) {
75 first_layer_clusters[i].add(first_layer_clusters[j]);
76 first_layer_clusters[j].clear();
77 merge = true;
78 }
79 }
80 // } else break;
81 }
82 } while (merge);
83 ldmx_log(trace) << "--- ELECTRON SEPARATION ---";
84 for (int i = 0; i < first_layer_clusters.size(); i++) {
85 if (first_layer_clusters[i].empty()) continue;
86 ldmx_log(trace) << " Cluster " << i
87 << " x: " << first_layer_clusters[i].centroidX()
88 << " y: " << first_layer_clusters[i].centroidY();
89 for (int j = i + 1; j < first_layer_clusters.size(); j++) {
90 if (first_layer_clusters[j].empty()) continue;
91 auto d = dist(first_layer_clusters[i].centroidX(),
92 first_layer_clusters[i].centroidY(),
93 first_layer_clusters[j].centroidX(),
94 first_layer_clusters[j].centroidY());
95 ldmx_log(trace) << "Dist to cluster " << j << ": " << d;
96 }
97 }
98}

◆ getCentroidDistances()

std::vector< double > ecal::CLUE::getCentroidDistances ( ) const
inline

Definition at line 102 of file CLUE.h.

102 {
103 return centroid_distances_;
104 }

◆ getClusters()

std::vector< IntermediateCluster > ecal::CLUE::getClusters ( ) const
inline

Definition at line 110 of file CLUE.h.

110 {
111 return final_clusters_;
112 }

◆ getFirstLayerCentroids()

std::vector< IntermediateCluster > ecal::CLUE::getFirstLayerCentroids ( ) const
inline

Definition at line 116 of file CLUE.h.

116 {
117 return first_layer_centroids_;
118 }

◆ getInitialClusterNbr()

int ecal::CLUE::getInitialClusterNbr ( ) const
inline

Definition at line 108 of file CLUE.h.

108{ return initial_cluster_nbr_; }

◆ getNLoops()

int ecal::CLUE::getNLoops ( ) const
inline

Definition at line 106 of file CLUE.h.

106{ return clustering_loops_; }

◆ roundToDecimal()

float ecal::CLUE::roundToDecimal ( float x,
int num_decimal_precision_digits )

Definition at line 153 of file CLUE.cxx.

153 {
154 float power_of_10 = std::pow(10, num_decimal_precision_digits);
155 return std::round(x_ * power_of_10) / power_of_10;
156}

◆ setup()

std::vector< std::shared_ptr< CLUE::Density > > ecal::CLUE::setup ( const std::vector< const ldmx::EcalHit * > & hits)

Definition at line 158 of file CLUE.cxx.

159 {
160 std::vector<std::shared_ptr<Density>> densities;
161 std::map<std::pair<float, float>, std::shared_ptr<Density>> density_map;
162 event_centroid_ = IntermediateCluster();
163 ldmx_log(trace) << "--- SETUP ---";
164 ldmx_log(trace) << "Building densities";
165 for (const auto& hit : hits) {
166 // collapse z dimension
167 float x = roundToDecimal(hit->getXPos(), 4);
168 float y = roundToDecimal(hit->getYPos(), 4);
169 float z = roundToDecimal(hit->getZPos(), 4);
170 ldmx_log(trace) << " New hit { x: " << x << " y: " << y << "}"
171 << " (and z: " << z << ")";
172 std::pair<float, float> coords;
173 if (dc_ != 0 && nbr_of_layers_ > 1) {
174 // if more than one layer, divide hit into densities with side dc
175 double i = std::ceil(std::abs(x) / dc_);
176 double j = std::ceil(std::abs(y) / dc_);
177 if (x < 0) {
178 i = -i;
179 x = (i + 0.5) * dc_;
180 } else {
181 x = (i - 0.5) * dc_;
182 }
183 if (y < 0) {
184 j = -j;
185 y = (j + 0.5) * dc_;
186 } else {
187 y = (j - 0.5) * dc_;
188 }
189 coords = {i, j};
190 ldmx_log(trace) << " Index " << i << ", " << j << "; x: " << x
191 << " y: " << y;
192 } else {
193 // if just one layer, have all densities with the same x,y be in same
194 // density
195 coords = {x, y};
196 }
197
198 if (density_map.find(coords) == density_map.end()) {
199 density_map.emplace(coords, std::make_shared<CLUE::Density>(x, y));
200 ldmx_log(trace) << " * New density created";
201 } else {
202 ldmx_log(trace) << " --> Found density with x: "
203 << density_map[coords]->x_
204 << " y: " << density_map[coords]->y_;
205 }
206 density_map[coords]->hits_.push_back(hit);
207 density_map[coords]->total_energy_ += hit->getEnergy();
208 density_map[coords]->z_ += hit->getZPos();
209
210 event_centroid_.add(hit);
211 } // end of loop over hits
212
213 densities.reserve(density_map.size());
214 for (const auto& entry : density_map) {
215 densities.push_back(std::move(entry.second));
216 }
217 // sort according to energy
218 std::sort(densities.begin(), densities.end(),
219 [](const std::shared_ptr<CLUE::Density>& a,
220 const std::shared_ptr<CLUE::Density>& b) {
221 return a->total_energy_ > b->total_energy_;
222 });
223
224 ldmx_log(trace) << "Decide parents";
225
226 // decide delta_ and follower_of_
227 for (int i = 0; i < densities.size(); i++) {
228 densities[i]->index_ = i;
229 // avg z position
230 densities[i]->z_ = densities[i]->z_ / densities[i]->hits_.size();
231 ldmx_log(trace) << " Index: " << i << "; x: " << densities[i]->x_
232 << "; y: " << densities[i]->y_
233 << "; Energy: " << densities[i]->total_energy_;
234 // loop through all higher energy densities
235 for (int j = 0; j < i; j++) {
236 float distance_2d = dist(densities[i]->x_, densities[i]->y_,
237 densities[j]->x_, densities[j]->y_);
238 // condition energyJ > energyI but this should be baked in as we sorted
239 // according to energy
240 if ((distance_2d < dm_) && (distance_2d < densities[i]->delta_)) {
241 densities[i]->delta_ = distance_2d;
242 densities[i]->follower_of_ = j;
243 ldmx_log(trace) << " New parent, index " << j
244 << "; delta_2d: " << std::setprecision(4)
245 << distance_2d;
246 }
247 }
248 }
249 return densities;
250} // end of setup
void add(const HitType *hit)
Add a hit to the cluster using its stored position.

◆ setupForClue3D()

std::vector< std::shared_ptr< CLUE::Density > > ecal::CLUE::setupForClue3D ( )

Definition at line 451 of file CLUE.cxx.

451 {
452 ldmx_log(trace) << "--- LAYER SETUP ---";
453 std::vector<std::shared_ptr<CLUE::Density>> densities;
454 layer_rho_c_.clear();
455 for (int layer = 0; layer < nbr_of_layers_; layer++) {
456 ldmx_log(trace) << " LAYER " << layer << " with " << seeds_[layer].size()
457 << " seeds";
458 auto& seeds_in_current_layer = seeds_[layer];
459 double highest_energy = 0.;
460 for (const auto& current_seed : seeds_in_current_layer) {
461 // for each seed in layer
462 current_seed->layer_ = layer;
463 if (current_seed->total_energy_ > highest_energy)
464 highest_energy = current_seed->total_energy_;
465 ldmx_log(trace) << " Density with index " << current_seed->index_
466 << ", energy: " << current_seed->total_energy_
467 << " position (x,y)= {" << current_seed->x_ << ","
468 << current_seed->y_ << ")";
469 int depth = 1;
470 // decide delta_ and followerof from seeds in previous and next layer_
471 // do {
472 // depth++;
473 if ((layer - depth >= 0) && (layer - depth < seeds_.size())) {
474 ldmx_log(trace) << " Looking at pre-layer: " << layer - depth;
475 // look at previous layer
476 ldmx_log(trace) << " In previous layer... ";
477 auto& previous_layer = seeds_[layer - depth];
478 ldmx_log(trace) << " Got " << previous_layer.size() << " seeds";
479 for (const auto& prev_seed : previous_layer) {
480 // for each seed in previous layer
481 auto distance_2d_prev = dist(current_seed->x_, current_seed->y_,
482 prev_seed->x_, prev_seed->y_);
483 auto dz_prev = std::abs(current_seed->z_ - prev_seed->z_);
484 ldmx_log(trace) << " DeltaXY to index " << prev_seed->index_
485 << ": " << std::setprecision(4) << distance_2d_prev;
486 ldmx_log(trace) << " DeltaZ to index " << prev_seed->index_
487 << ": " << std::setprecision(4) << dz_prev;
488 if (prev_seed->total_energy_ > current_seed->total_energy_ &&
489 distance_2d_prev < current_seed->delta_ &&
490 dz_prev < current_seed->z_delta_) {
491 ldmx_log(trace) << " New parent: index " << prev_seed->index_
492 << " on layer " << layer - depth << "; energy "
493 << prev_seed->total_energy_;
494 ldmx_log(trace) << " New 2D distance to prev-layer: "
495 << std::setprecision(4) << distance_2d_prev;
496 ldmx_log(trace)
497 << " New delta Z: " << std::setprecision(4) << dz_prev;
498 current_seed->delta_ = distance_2d_prev;
499 current_seed->z_delta_ = dz_prev;
500 current_seed->follower_of_ = prev_seed->index_;
501 } else if (prev_seed->total_energy_ < current_seed->total_energy_) {
502 ldmx_log(trace) << " Breaking on index " << prev_seed->index_
503 << " with energy " << prev_seed->total_energy_;
504 break;
505 }
506 }
507 }
508
509 if (layer + depth < nbr_of_layers_ && layer + depth < seeds_.size()) {
510 ldmx_log(trace) << " Looking at post-layer: " << layer + depth;
511 auto& next_layer = seeds_[layer + depth];
512 ldmx_log(trace) << " Got " << next_layer.size() << " seeds";
513 for (const auto& next_seed : next_layer) {
514 auto distance_2d_next = dist(current_seed->x_, current_seed->y_,
515 next_seed->x_, next_seed->y_);
516 auto dz_next = std::abs(current_seed->z_ - next_seed->z_);
517 ldmx_log(trace) << " DeltaXY to index " << next_seed->index_
518 << ": " << std::setprecision(4) << distance_2d_next;
519 ldmx_log(trace) << " DeltaZ to index_" << next_seed->index_
520 << ": " << std::setprecision(4) << dz_next;
521 if (next_seed->total_energy_ > current_seed->total_energy_ &&
522 distance_2d_next < current_seed->delta_ &&
523 dz_next < current_seed->z_delta_) {
524 ldmx_log(trace) << " New parent: index_" << next_seed->index_
525 << " on layer " << layer + depth << "; energy "
526 << next_seed->total_energy_;
527 ldmx_log(trace) << " New 2D distance to next-layer: "
528 << std::setprecision(4) << distance_2d_next;
529 ldmx_log(trace)
530 << " New delta_Z: " << std::setprecision(4) << dz_next;
531 current_seed->delta_ = distance_2d_next;
532 current_seed->z_delta_ = dz_next;
533 current_seed->follower_of_ = next_seed->index_;
534 } else if (next_seed->total_energy_ < current_seed->total_energy_) {
535 ldmx_log(trace) << " Breaking on index_" << next_seed->index_
536 << " with energy " << next_seed->total_energy_;
537 break;
538 }
539 } // end loop on seeds in next layer
540 } // end of looking at next layer
541 // } while (currentLayer[i]->layerFollowerOf == -1 && (layer - depth >=
542 // 0 || layer + depth < nbr_of_layers_));
543 ldmx_log(trace) << " Done setting parents.";
544 densities.push_back(current_seed);
545 }
546 // TODO: This 2 needs to be configurable
547 layer_rho_c_.push_back(highest_energy / 2);
548 } // end loop over layers
549 return densities;
550} // end of setupForClue3D

Member Data Documentation

◆ centroid_distances_

std::vector<double> ecal::CLUE::centroid_distances_
private

Definition at line 150 of file CLUE.h.

◆ clustering_loops_

int ecal::CLUE::clustering_loops_
private

Definition at line 121 of file CLUE.h.

◆ dc_

double ecal::CLUE::dc_
private

Definition at line 125 of file CLUE.h.

◆ deltac_

double ecal::CLUE::deltac_
private

Definition at line 127 of file CLUE.h.

◆ deltao_

double ecal::CLUE::deltao_
private

Definition at line 128 of file CLUE.h.

◆ dm_

double ecal::CLUE::dm_
private

Definition at line 129 of file CLUE.h.

◆ event_centroid_

IntermediateCluster ecal::CLUE::event_centroid_
private

Definition at line 151 of file CLUE.h.

◆ final_clusters_

std::vector<IntermediateCluster> ecal::CLUE::final_clusters_
private

Definition at line 159 of file CLUE.h.

◆ first_layer_centroids_

std::vector<IntermediateCluster> ecal::CLUE::first_layer_centroids_
private

Definition at line 153 of file CLUE.h.

◆ initial_cluster_nbr_

int ecal::CLUE::initial_cluster_nbr_ {-1}
private

Definition at line 158 of file CLUE.h.

158{-1};

◆ layer_centroid_separations_

std::vector<std::pair<double, double> > ecal::CLUE::layer_centroid_separations_
private

Definition at line 160 of file CLUE.h.

◆ layer_delta_c_

std::vector<double> ecal::CLUE::layer_delta_c_
private

Definition at line 136 of file CLUE.h.

◆ layer_rho_c_

std::vector<double> ecal::CLUE::layer_rho_c_
private

Definition at line 135 of file CLUE.h.

◆ max_layers_

int ecal::CLUE::max_layers_ {32}
private

Definition at line 132 of file CLUE.h.

132{32};

◆ nbr_of_layers_

int ecal::CLUE::nbr_of_layers_
private

Definition at line 133 of file CLUE.h.

◆ radius_

std::vector<double> ecal::CLUE::radius_
private
Initial value:
{
5.723387467629167, 5.190678018534044, 5.927290663506518,
6.182560329200212, 7.907549398117859, 8.606100542857211,
10.93381822596916, 12.043201938160239, 14.784548371508041,
16.102403056546482, 18.986402399412817, 20.224453740305716,
23.048820910305643, 24.11202594672678, 26.765135236851666,
27.78700483852502, 30.291794353801293, 31.409870873194464,
33.91006482486666, 35.173073672355926, 38.172422630271,
40.880288341493205, 44.696485719120005, 49.23802839743545,
53.789910813378675, 60.87843355562641, 66.32931132415688,
75.78117972604727, 86.04697356716805, 96.90360704034346}

Definition at line 138 of file CLUE.h.

138 {
139 5.723387467629167, 5.190678018534044, 5.927290663506518,
140 6.182560329200212, 7.907549398117859, 8.606100542857211,
141 10.93381822596916, 12.043201938160239, 14.784548371508041,
142 16.102403056546482, 18.986402399412817, 20.224453740305716,
143 23.048820910305643, 24.11202594672678, 26.765135236851666,
144 27.78700483852502, 30.291794353801293, 31.409870873194464,
145 33.91006482486666, 35.173073672355926, 38.172422630271,
146 40.880288341493205, 44.696485719120005, 49.23802839743545,
147 53.789910813378675, 60.87843355562641, 66.32931132415688,
148 75.78117972604727, 86.04697356716805, 96.90360704034346};

◆ reclustering_

bool ecal::CLUE::reclustering_
private

Definition at line 123 of file CLUE.h.

◆ rhoc_

double ecal::CLUE::rhoc_
private

Definition at line 126 of file CLUE.h.

◆ seed_index_

int ecal::CLUE::seed_index_ {0}
private

Definition at line 155 of file CLUE.h.

155{0};

◆ seeds_

std::vector<std::vector<std::shared_ptr<Density> > > ecal::CLUE::seeds_
private

Definition at line 156 of file CLUE.h.


The documentation for this class was generated from the following files: