CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
LCToSCAssociatorByEnergyScoreImpl Class Reference

#include <LCToSCAssociatorByEnergyScoreImpl.h>

Inheritance diagram for LCToSCAssociatorByEnergyScoreImpl:
hgcal::LayerClusterToSimClusterAssociatorBaseImpl

Public Member Functions

hgcal::RecoToSimCollectionWithSimClusters associateRecoToSim (const edm::Handle< reco::CaloClusterCollection > &cCH, const edm::Handle< SimClusterCollection > &sCCH) const override
 Associate a LayerCluster to SimClusters. More...
 
hgcal::SimToRecoCollectionWithSimClusters associateSimToReco (const edm::Handle< reco::CaloClusterCollection > &cCH, const edm::Handle< SimClusterCollection > &sCCH) const override
 Associate a SimCluster to LayerClusters. More...
 
 LCToSCAssociatorByEnergyScoreImpl (edm::EDProductGetter const &, bool, std::shared_ptr< hgcal::RecHitTools >, const std::unordered_map< DetId, const HGCRecHit *> *)
 
- Public Member Functions inherited from hgcal::LayerClusterToSimClusterAssociatorBaseImpl
 LayerClusterToSimClusterAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~LayerClusterToSimClusterAssociatorBaseImpl ()
 Destructor. More...
 

Private Member Functions

hgcal::association makeConnections (const edm::Handle< reco::CaloClusterCollection > &cCH, const edm::Handle< SimClusterCollection > &sCCH) const
 

Private Attributes

const bool hardScatterOnly_
 
const std::unordered_map< DetId, const HGCRecHit * > * hitMap_
 
unsigned layers_
 
edm::EDProductGetter const * productGetter_
 
std::shared_ptr< hgcal::RecHitToolsrecHitTools_
 

Detailed Description

Definition at line 62 of file LCToSCAssociatorByEnergyScoreImpl.h.

Constructor & Destructor Documentation

◆ LCToSCAssociatorByEnergyScoreImpl()

LCToSCAssociatorByEnergyScoreImpl::LCToSCAssociatorByEnergyScoreImpl ( edm::EDProductGetter const &  productGetter,
bool  hardScatterOnly,
std::shared_ptr< hgcal::RecHitTools recHitTools,
const std::unordered_map< DetId, const HGCRecHit *> *  hitMap 
)
explicit

Definition at line 6 of file LCToSCAssociatorByEnergyScoreImpl.cc.

References layers_, and recHitTools_.

11  : hardScatterOnly_(hardScatterOnly), recHitTools_(recHitTools), hitMap_(hitMap), productGetter_(&productGetter) {
12  layers_ = recHitTools_->lastLayerBH();
13 }
const std::unordered_map< DetId, const HGCRecHit * > * hitMap_
std::shared_ptr< hgcal::RecHitTools > recHitTools_
EDProductGetter const * productGetter(std::atomic< void const *> const &iCache)

Member Function Documentation

◆ associateRecoToSim()

hgcal::RecoToSimCollectionWithSimClusters LCToSCAssociatorByEnergyScoreImpl::associateRecoToSim ( const edm::Handle< reco::CaloClusterCollection > &  cCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
overridevirtual

Associate a LayerCluster to SimClusters.

Reimplemented from hgcal::LayerClusterToSimClusterAssociatorBaseImpl.

Definition at line 499 of file LCToSCAssociatorByEnergyScoreImpl.cc.

References edm::AssociationMap< Tag >::insert(), electronStore::links, LogDebug, makeConnections(), and productGetter_.

500  {
502  const auto& links = makeConnections(cCCH, sCCH);
503 
504  const auto& scsInLayerCluster = std::get<0>(links);
505  for (size_t lcId = 0; lcId < scsInLayerCluster.size(); ++lcId) {
506  for (auto& scPair : scsInLayerCluster[lcId]) {
507  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
508  << "layerCluster Id: \t" << lcId << "\t SC id: \t" << scPair.first << "\t score \t" << scPair.second << "\n";
509  // Fill AssociationMap
510  returnValue.insert(edm::Ref<reco::CaloClusterCollection>(cCCH, lcId), // Ref to LC
511  std::make_pair(edm::Ref<SimClusterCollection>(sCCH, scPair.first),
512  scPair.second) // Pair <Ref to SC, score>
513  );
514  }
515  }
516  return returnValue;
517 }
hgcal::association makeConnections(const edm::Handle< reco::CaloClusterCollection > &cCH, const edm::Handle< SimClusterCollection > &sCCH) const
#define LogDebug(id)

◆ associateSimToReco()

hgcal::SimToRecoCollectionWithSimClusters LCToSCAssociatorByEnergyScoreImpl::associateSimToReco ( const edm::Handle< reco::CaloClusterCollection > &  cCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
overridevirtual

Associate a SimCluster to LayerClusters.

Reimplemented from hgcal::LayerClusterToSimClusterAssociatorBaseImpl.

Definition at line 519 of file LCToSCAssociatorByEnergyScoreImpl.cc.

References edm::AssociationMap< Tag >::insert(), electronStore::links, makeConnections(), and productGetter_.

520  {
522  const auto& links = makeConnections(cCCH, sCCH);
523  const auto& lcsInSimCluster = std::get<1>(links);
524  for (size_t scId = 0; scId < lcsInSimCluster.size(); ++scId) {
525  for (size_t layerId = 0; layerId < lcsInSimCluster[scId].size(); ++layerId) {
526  for (auto& lcPair : lcsInSimCluster[scId][layerId].layerClusterIdToEnergyAndScore) {
527  returnValue.insert(
528  edm::Ref<SimClusterCollection>(sCCH, scId), // Ref to SC
529  std::make_pair(edm::Ref<reco::CaloClusterCollection>(cCCH, lcPair.first), // Pair <Ref to LC,
530  std::make_pair(lcPair.second.first, lcPair.second.second)) // pair <energy, score> >
531  );
532  }
533  }
534  }
535  return returnValue;
536 }
hgcal::association makeConnections(const edm::Handle< reco::CaloClusterCollection > &cCH, const edm::Handle< SimClusterCollection > &sCCH) const

◆ makeConnections()

hgcal::association LCToSCAssociatorByEnergyScoreImpl::makeConnections ( const edm::Handle< reco::CaloClusterCollection > &  cCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
private

Definition at line 15 of file LCToSCAssociatorByEnergyScoreImpl.cc.

References c, bsc_activity_cfg::clusters, relativeConstraints::empty, mps_fire::end, HCALHighEnergyHPDFilter_cfi::energy, edmPickEvents::event, f, spr::find(), cms::cuda::for(), h, hardScatterOnly_, hitMap_, mps_fire::i, dqmiolumiharvest::j, dqmdumpme::last, layers_, LogDebug, SiStripPI::min, or, funct::pow(), edm::Handle< T >::product(), recHitTools_, jetUpdater_cfi::sort, and tier0::unique().

Referenced by associateRecoToSim(), and associateSimToReco().

16  {
17  // Get collections
18  const auto& clusters = *cCCH.product();
19  const auto& simClusters = *sCCH.product();
20  auto nLayerClusters = clusters.size();
21 
22  //There shouldn't be any SimTracks from different crossings, but maybe they will be added later.
23  //At the moment there should be one SimTrack in each SimCluster.
24  auto nSimClusters = simClusters.size();
25  std::vector<size_t> sCIndices;
26  for (unsigned int scId = 0; scId < nSimClusters; ++scId) {
27  if (hardScatterOnly_ && (simClusters[scId].g4Tracks()[0].eventId().event() != 0 or
28  simClusters[scId].g4Tracks()[0].eventId().bunchCrossing() != 0)) {
29  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
30  << "Excluding SimCluster from event: " << simClusters[scId].g4Tracks()[0].eventId().event()
31  << " with BX: " << simClusters[scId].g4Tracks()[0].eventId().bunchCrossing() << std::endl;
32  continue;
33  }
34  sCIndices.emplace_back(scId);
35  }
36  nSimClusters = sCIndices.size();
37 
38  // Initialize lcsInSimCluster. It contains the simClusterOnLayer structure for all simClusters in each layer and
39  // among other the information to compute the SimCluster-To-LayerCluster score. It is one of the two objects that
40  // build the output of the makeConnections function.
41  // lcsInSimCluster[scId][layerId]
42  hgcal::simClusterToLayerCluster lcsInSimCluster;
43  lcsInSimCluster.resize(nSimClusters);
44  for (unsigned int i = 0; i < nSimClusters; ++i) {
45  lcsInSimCluster[i].resize(layers_ * 2);
46  for (unsigned int j = 0; j < layers_ * 2; ++j) {
47  lcsInSimCluster[i][j].simClusterId = i;
48  lcsInSimCluster[i][j].energy = 0.f;
49  lcsInSimCluster[i][j].hits_and_fractions.clear();
50  }
51  }
52 
53  // Fill detIdToSimClusterId_Map and update lcsInSimCluster
54  // The detIdToSimClusterId_Map is used to connect a hit Detid (key) with all the SimClusters that
55  // contributed to that hit by storing the SimCluster id and the fraction of the hit. Observe here
56  // that in contrast to the CaloParticle case there is no merging and summing of the fractions, which
57  // in the CaloParticle's case was necessary due to the multiple SimClusters of a single CaloParticle.
58  std::unordered_map<DetId, std::vector<hgcal::detIdInfoInCluster>> detIdToSimClusterId_Map;
59  for (const auto& scId : sCIndices) {
60  const auto& hits_and_fractions = simClusters[scId].hits_and_fractions();
61  for (const auto& it_haf : hits_and_fractions) {
62  const auto hitid = (it_haf.first);
63  const auto scLayerId =
64  recHitTools_->getLayerWithOffset(hitid) + layers_ * ((recHitTools_->zside(hitid) + 1) >> 1) - 1;
65 
66  const auto itcheck = hitMap_->find(hitid);
67  if (itcheck != hitMap_->end()) {
68  auto hit_find_it = detIdToSimClusterId_Map.find(hitid);
69  if (hit_find_it == detIdToSimClusterId_Map.end()) {
70  detIdToSimClusterId_Map[hitid] = std::vector<hgcal::detIdInfoInCluster>();
71  }
72  detIdToSimClusterId_Map[hitid].emplace_back(scId, it_haf.second);
73 
74  const HGCRecHit* hit = itcheck->second;
75  lcsInSimCluster[scId][scLayerId].energy += it_haf.second * hit->energy();
76  lcsInSimCluster[scId][scLayerId].hits_and_fractions.emplace_back(hitid, it_haf.second);
77  }
78  }
79  } // end of loop over SimClusters
80 
81 #ifdef EDM_ML_DEBUG
82  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
83  << "lcsInSimCluster INFO (Only SimCluster filled at the moment)" << std::endl;
84  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << " # of clusters : " << nLayerClusters << std::endl;
85  for (size_t sc = 0; sc < lcsInSimCluster.size(); ++sc) {
86  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << "For SimCluster Idx: " << sc << " we have: " << std::endl;
87  for (size_t sclay = 0; sclay < lcsInSimCluster[sc].size(); ++sclay) {
88  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << " On Layer: " << sclay << " we have:" << std::endl;
89  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
90  << " SimClusterIdx: " << lcsInSimCluster[sc][sclay].simClusterId << std::endl;
91  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
92  << " Energy: " << lcsInSimCluster[sc][sclay].energy << std::endl;
93  double tot_energy = 0.;
94  for (auto const& haf : lcsInSimCluster[sc][sclay].hits_and_fractions) {
95  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
96  << " Hits/fraction/energy: " << (uint32_t)haf.first << "/" << haf.second << "/"
97  << haf.second * hitMap_->at(haf.first)->energy() << std::endl;
98  tot_energy += haf.second * hitMap_->at(haf.first)->energy();
99  }
100  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << " Tot Sum haf: " << tot_energy << std::endl;
101  for (auto const& lc : lcsInSimCluster[sc][sclay].layerClusterIdToEnergyAndScore) {
102  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << " lcIdx/energy/score: " << lc.first << "/"
103  << lc.second.first << "/" << lc.second.second << std::endl;
104  }
105  }
106  }
107 
108  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << "detIdToSimClusterId_Map INFO" << std::endl;
109  for (auto const& sc : detIdToSimClusterId_Map) {
110  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
111  << "For detId: " << (uint32_t)sc.first
112  << " we have found the following connections with SimClusters:" << std::endl;
113  // At this point here if you activate the printing you will notice cases where in a
114  // specific detId there are more that one SimClusters contributing with fractions less than 1.
115  // This is important since it effects the score computation, since the fraction is also in the
116  // denominator of the score formula.
117  for (auto const& sclu : sc.second) {
118  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
119  << " SimCluster Id: " << sclu.clusterId << " with fraction: " << sclu.fraction
120  << " and energy: " << sclu.fraction * hitMap_->at(sc.first)->energy() << std::endl;
121  }
122  }
123 #endif
124 
125  // Fill detIdToLayerClusterId_Map and scsInLayerCluster; update lcsInSimCluster
126  // The detIdToLayerClusterId_Map is used to connect a hit Detid (key) with all the LayerClusters that
127  // contributed to that hit by storing the LayerCluster id and the fraction of the corresponding hit.
128  std::unordered_map<DetId, std::vector<hgcal::detIdInfoInCluster>> detIdToLayerClusterId_Map;
129  // scsInLayerCluster together with lcsInSimCluster are the two objects that are used to build the
130  // output of the makeConnections function. scsInLayerCluster connects a LayerCluster with
131  // all the SimClusters that share at least one cell with the LayerCluster and for each pair (LC,SC)
132  // it stores the score.
133  hgcal::layerClusterToSimCluster scsInLayerCluster; //[lcId][scId]->(score)
134  scsInLayerCluster.resize(nLayerClusters);
135 
136  for (unsigned int lcId = 0; lcId < nLayerClusters; ++lcId) {
137  const std::vector<std::pair<DetId, float>>& hits_and_fractions = clusters[lcId].hitsAndFractions();
138  unsigned int numberOfHitsInLC = hits_and_fractions.size();
139  const auto firstHitDetId = hits_and_fractions[0].first;
140  int lcLayerId =
141  recHitTools_->getLayerWithOffset(firstHitDetId) + layers_ * ((recHitTools_->zside(firstHitDetId) + 1) >> 1) - 1;
142 
143  for (unsigned int hitId = 0; hitId < numberOfHitsInLC; hitId++) {
144  const auto rh_detid = hits_and_fractions[hitId].first;
145  const auto rhFraction = hits_and_fractions[hitId].second;
146 
147  auto hit_find_in_LC = detIdToLayerClusterId_Map.find(rh_detid);
148  if (hit_find_in_LC == detIdToLayerClusterId_Map.end()) {
149  detIdToLayerClusterId_Map[rh_detid] = std::vector<hgcal::detIdInfoInCluster>();
150  }
151  detIdToLayerClusterId_Map[rh_detid].emplace_back(lcId, rhFraction);
152 
153  auto hit_find_in_SC = detIdToSimClusterId_Map.find(rh_detid);
154 
155  if (hit_find_in_SC != detIdToSimClusterId_Map.end()) {
156  const auto itcheck = hitMap_->find(rh_detid);
157  const HGCRecHit* hit = itcheck->second;
158  //Loops through all the simclusters that have the layer cluster rechit under study
159  //Here is time to update the lcsInSimCluster and connect the SimCluster with all
160  //the layer clusters that have the current rechit detid matched.
161  for (auto& h : hit_find_in_SC->second) {
162  //lcsInSimCluster[simclusterId][layerId][layerclusterId]-> (energy,score)
163  //SC_i - > LC_j, LC_k, ...
164  lcsInSimCluster[h.clusterId][lcLayerId].layerClusterIdToEnergyAndScore[lcId].first +=
165  h.fraction * hit->energy();
166  //LC_i -> SC_j, SC_k, ...
167  scsInLayerCluster[lcId].emplace_back(h.clusterId, 0.f);
168  }
169  }
170  } // End loop over hits on a LayerCluster
171  } // End of loop over LayerClusters
172 
173 #ifdef EDM_ML_DEBUG
174  for (unsigned int lcId = 0; lcId < nLayerClusters; ++lcId) {
175  const auto& hits_and_fractions = clusters[lcId].hitsAndFractions();
176  unsigned int numberOfHitsInLC = hits_and_fractions.size();
177  const auto firstHitDetId = hits_and_fractions[0].first;
178  int lcLayerId =
179  recHitTools_->getLayerWithOffset(firstHitDetId) + layers_ * ((recHitTools_->zside(firstHitDetId) + 1) >> 1) - 1;
180 
181  // This vector will store, for each hit in the Layercluster, the index of
182  // the SimCluster that contributed the most, in terms of energy, to it.
183  // Special values are:
184  //
185  // -2 --> the reconstruction fraction of the RecHit is 0 (used in the past to monitor Halo Hits)
186  // -3 --> same as before with the added condition that no SimCluster has been linked to this RecHit
187  // -1 --> the reco fraction is >0, but no SimCluster has been linked to it
188  // >=0 --> index of the linked SimCluster
189  std::vector<int> hitsToSimClusterId(numberOfHitsInLC);
190  // This will store the index of the SimCluster linked to the LayerCluster that has the most number of hits in common.
191  int maxSCId_byNumberOfHits = -1;
192  // This will store the maximum number of shared hits between a Layercluster and a SimCluster
193  unsigned int maxSCNumberOfHitsInLC = 0;
194  // This will store the index of the SimCluster linked to the LayerCluster that has the most energy in common.
195  int maxSCId_byEnergy = -1;
196  // This will store the maximum number of shared energy between a Layercluster and a SimCluster
197  float maxEnergySharedLCandSC = 0.f;
198  // This will store the fraction of the LayerCluster energy shared with the best(energy) SimCluster: e_shared/lc_energy
199  float energyFractionOfLCinSC = 0.f;
200  // This will store the fraction of the SimCluster energy shared with the LayerCluster: e_shared/sc_energy
201  float energyFractionOfSCinLC = 0.f;
202  std::unordered_map<unsigned, unsigned> occurrencesSCinLC;
203  unsigned int numberOfNoiseHitsInLC = 0;
204  std::unordered_map<unsigned, float> SCEnergyInLC;
205 
206  for (unsigned int hitId = 0; hitId < numberOfHitsInLC; hitId++) {
207  const auto rh_detid = hits_and_fractions[hitId].first;
208  const auto rhFraction = hits_and_fractions[hitId].second;
209 
210  auto hit_find_in_SC = detIdToSimClusterId_Map.find(rh_detid);
211 
212  // if the fraction is zero or the hit does not belong to any simcluster
213  // set the caloparticleId for the hit to -1 this will
214  // contribute to the number of noise hits
215 
216  // MR Remove the case in which the fraction is 0, since this could be a
217  // real hit that has been marked as halo.
218  if (rhFraction == 0.) {
219  hitsToSimClusterId[hitId] = -2;
220  }
221  //Now check if there are SimClusters linked to this rechit of the layercluster
222  if (hit_find_in_SC == detIdToSimClusterId_Map.end()) {
223  hitsToSimClusterId[hitId] -= 1;
224  } else {
225  const auto itcheck = hitMap_->find(rh_detid);
226  const HGCRecHit* hit = itcheck->second;
227  auto maxSCEnergyInLC = 0.f;
228  auto maxSCId = -1;
229  //Loop through all the linked SimClusters
230  for (auto& h : hit_find_in_SC->second) {
231  SCEnergyInLC[h.clusterId] += h.fraction * hit->energy();
232  // Keep track of which SimCluster ccontributed the most, in terms
233  // of energy, to this specific LayerCluster.
234  if (SCEnergyInLC[h.clusterId] > maxSCEnergyInLC) {
235  maxSCEnergyInLC = SCEnergyInLC[h.clusterId];
236  maxSCId = h.clusterId;
237  }
238  }
239  hitsToSimClusterId[hitId] = maxSCId;
240  }
241  } // End loop over hits on a LayerCluster
242 
243  for (const auto& c : hitsToSimClusterId) {
244  if (c < 0) {
245  numberOfNoiseHitsInLC++;
246  } else {
247  occurrencesSCinLC[c]++;
248  }
249  }
250 
251  for (const auto& c : occurrencesSCinLC) {
252  if (c.second > maxSCNumberOfHitsInLC) {
253  maxSCId_byNumberOfHits = c.first;
254  maxSCNumberOfHitsInLC = c.second;
255  }
256  }
257 
258  for (const auto& c : SCEnergyInLC) {
259  if (c.second > maxEnergySharedLCandSC) {
260  maxSCId_byEnergy = c.first;
261  maxEnergySharedLCandSC = c.second;
262  }
263  }
264 
265  float totalSCEnergyOnLayer = 0.f;
266  if (maxSCId_byEnergy >= 0) {
267  totalSCEnergyOnLayer = lcsInSimCluster[maxSCId_byEnergy][lcLayerId].energy;
268  energyFractionOfSCinLC = maxEnergySharedLCandSC / totalSCEnergyOnLayer;
269  if (clusters[lcId].energy() > 0.f) {
270  energyFractionOfLCinSC = maxEnergySharedLCandSC / clusters[lcId].energy();
271  }
272  }
273 
274  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << std::setw(10) << "LayerId:"
275  << "\t" << std::setw(12) << "layerCluster"
276  << "\t" << std::setw(10) << "lc energy"
277  << "\t" << std::setw(5) << "nhits"
278  << "\t" << std::setw(12) << "noise hits"
279  << "\t" << std::setw(22) << "maxSCId_byNumberOfHits"
280  << "\t" << std::setw(8) << "nhitsSC"
281  << "\t" << std::setw(13) << "maxSCId_byEnergy"
282  << "\t" << std::setw(20) << "maxEnergySharedLCandSC"
283  << "\t" << std::setw(22) << "totalSCEnergyOnLayer"
284  << "\t" << std::setw(22) << "energyFractionOfLCinSC"
285  << "\t" << std::setw(25) << "energyFractionOfSCinLC"
286  << "\t"
287  << "\n";
288  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
289  << std::setw(10) << lcLayerId << "\t" << std::setw(12) << lcId << "\t" << std::setw(10)
290  << clusters[lcId].energy() << "\t" << std::setw(5) << numberOfHitsInLC << "\t" << std::setw(12)
291  << numberOfNoiseHitsInLC << "\t" << std::setw(22) << maxSCId_byNumberOfHits << "\t" << std::setw(8)
292  << maxSCNumberOfHitsInLC << "\t" << std::setw(13) << maxSCId_byEnergy << "\t" << std::setw(20)
293  << maxEnergySharedLCandSC << "\t" << std::setw(22) << totalSCEnergyOnLayer << "\t" << std::setw(22)
294  << energyFractionOfLCinSC << "\t" << std::setw(25) << energyFractionOfSCinLC << "\n";
295  } // End of loop over LayerClusters
296 
297  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
298  << "Improved lcsInSimCluster INFO (Now containing the linked layer clusters id and energy - score still empty)"
299  << std::endl;
300  for (size_t sc = 0; sc < lcsInSimCluster.size(); ++sc) {
301  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << "For SimCluster Idx: " << sc << " we have: " << std::endl;
302  for (size_t sclay = 0; sclay < lcsInSimCluster[sc].size(); ++sclay) {
303  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << " On Layer: " << sclay << " we have:" << std::endl;
304  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
305  << " SimClusterIdx: " << lcsInSimCluster[sc][sclay].simClusterId << std::endl;
306  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
307  << " Energy: " << lcsInSimCluster[sc][sclay].energy << std::endl;
308  double tot_energy = 0.;
309  for (auto const& haf : lcsInSimCluster[sc][sclay].hits_and_fractions) {
310  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
311  << " Hits/fraction/energy: " << (uint32_t)haf.first << "/" << haf.second << "/"
312  << haf.second * hitMap_->at(haf.first)->energy() << std::endl;
313  tot_energy += haf.second * hitMap_->at(haf.first)->energy();
314  }
315  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << " Tot Sum haf: " << tot_energy << std::endl;
316  for (auto const& lc : lcsInSimCluster[sc][sclay].layerClusterIdToEnergyAndScore) {
317  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << " lcIdx/energy/score: " << lc.first << "/"
318  << lc.second.first << "/" << lc.second.second << std::endl;
319  }
320  }
321  }
322 
323  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << "Improved detIdToSimClusterId_Map INFO" << std::endl;
324  for (auto const& sc : detIdToSimClusterId_Map) {
325  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
326  << "For detId: " << (uint32_t)sc.first
327  << " we have found the following connections with SimClusters:" << std::endl;
328  for (auto const& sclu : sc.second) {
329  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
330  << " SimCluster Id: " << sclu.clusterId << " with fraction: " << sclu.fraction
331  << " and energy: " << sclu.fraction * hitMap_->at(sc.first)->energy() << std::endl;
332  }
333  }
334 #endif
335 
336  // Update scsInLayerCluster; compute the score LayerCluster-to-SimCluster,
337  // together with the returned AssociationMap
338  for (unsigned int lcId = 0; lcId < nLayerClusters; ++lcId) {
339  // The simclusters contributing to the layer clusters should already be unique.
340  // find the unique simclusters id contributing to the layer clusters
341  std::sort(scsInLayerCluster[lcId].begin(), scsInLayerCluster[lcId].end());
342  auto last = std::unique(scsInLayerCluster[lcId].begin(), scsInLayerCluster[lcId].end());
343  scsInLayerCluster[lcId].erase(last, scsInLayerCluster[lcId].end());
344  const auto& hits_and_fractions = clusters[lcId].hitsAndFractions();
345  unsigned int numberOfHitsInLC = hits_and_fractions.size();
346  // If a reconstructed LayerCluster has energy 0 but is linked to a
347  // SimCluster, assigned score 1
348  if (clusters[lcId].energy() == 0. && !scsInLayerCluster[lcId].empty()) {
349  for (auto& scPair : scsInLayerCluster[lcId]) {
350  scPair.second = 1.;
351  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << "layerClusterId : \t " << lcId << "\t SC id : \t"
352  << scPair.first << "\t score \t " << scPair.second << "\n";
353  }
354  continue;
355  }
356 
357  // Compute the correct normalization.
358  // It is the inverse of the denominator of the LCToSC score formula. Observe that this is the sum of the squares.
359  float invLayerClusterEnergyWeight = 0.f;
360  for (auto const& haf : hits_and_fractions) {
361  invLayerClusterEnergyWeight +=
362  (haf.second * hitMap_->at(haf.first)->energy()) * (haf.second * hitMap_->at(haf.first)->energy());
363  }
364  invLayerClusterEnergyWeight = 1.f / invLayerClusterEnergyWeight;
365  for (unsigned int i = 0; i < numberOfHitsInLC; ++i) {
366  DetId rh_detid = hits_and_fractions[i].first;
367  float rhFraction = hits_and_fractions[i].second;
368 
369  bool hitWithSC = (detIdToSimClusterId_Map.find(rh_detid) != detIdToSimClusterId_Map.end());
370 
371  auto itcheck = hitMap_->find(rh_detid);
372  const HGCRecHit* hit = itcheck->second;
373  float hitEnergyWeight = hit->energy() * hit->energy();
374 
375  for (auto& scPair : scsInLayerCluster[lcId]) {
376  float scFraction = 0.f;
377  if (hitWithSC) {
378  auto findHitIt = std::find(detIdToSimClusterId_Map[rh_detid].begin(),
379  detIdToSimClusterId_Map[rh_detid].end(),
380  hgcal::detIdInfoInCluster{scPair.first, 0.f});
381  if (findHitIt != detIdToSimClusterId_Map[rh_detid].end())
382  scFraction = findHitIt->fraction;
383  }
384  scPair.second += std::min(std::pow(rhFraction - scFraction, 2), std::pow(rhFraction, 2)) * hitEnergyWeight *
385  invLayerClusterEnergyWeight;
386 #ifdef EDM_ML_DEBUG
387  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
388  << "rh_detid:\t" << (uint32_t)rh_detid << "\tlayerClusterId:\t" << lcId << "\t"
389  << "rhfraction,scfraction:\t" << rhFraction << ", " << scFraction << "\t"
390  << "hitEnergyWeight:\t" << hitEnergyWeight << "\t"
391  << "current score:\t" << scPair.second << "\t"
392  << "invLayerClusterEnergyWeight:\t" << invLayerClusterEnergyWeight << "\n";
393 #endif
394  }
395  } // End of loop over Hits within a LayerCluster
396 #ifdef EDM_ML_DEBUG
397  if (scsInLayerCluster[lcId].empty())
398  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << "layerCluster Id: \t" << lcId << "\tSC id:\t-1 "
399  << "\t score \t-1"
400  << "\n";
401 #endif
402  } // End of loop over LayerClusters
403 
404  // Compute the SimCluster-To-LayerCluster score
405  for (const auto& scId : sCIndices) {
406  for (unsigned int layerId = 0; layerId < layers_ * 2; ++layerId) {
407  unsigned int SCNumberOfHits = lcsInSimCluster[scId][layerId].hits_and_fractions.size();
408  if (SCNumberOfHits == 0)
409  continue;
410 #ifdef EDM_ML_DEBUG
411  int lcWithMaxEnergyInSC = -1;
412  //energy of the most energetic LC from all that were linked to SC
413  float maxEnergyLCinSC = 0.f;
414  //Energy of the SC scId on layer layerId that was reconstructed.
415  float SCenergy = lcsInSimCluster[scId][layerId].energy;
416  //most energetic LC from all LCs linked to SC over SC energy.
417  float SCEnergyFractionInLC = 0.f;
418  for (auto& lc : lcsInSimCluster[scId][layerId].layerClusterIdToEnergyAndScore) {
419  if (lc.second.first > maxEnergyLCinSC) {
420  maxEnergyLCinSC = lc.second.first;
421  lcWithMaxEnergyInSC = lc.first;
422  }
423  }
424  if (SCenergy > 0.f)
425  SCEnergyFractionInLC = maxEnergyLCinSC / SCenergy;
426 
427  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
428  << std::setw(8) << "LayerId:\t" << std::setw(12) << "simcluster\t" << std::setw(15) << "sc total energy\t"
429  << std::setw(15) << "scEnergyOnLayer\t" << std::setw(14) << "SCNhitsOnLayer\t" << std::setw(18)
430  << "lcWithMaxEnergyInSC\t" << std::setw(15) << "maxEnergyLCinSC\t" << std::setw(20) << "SCEnergyFractionInLC"
431  << "\n";
432  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
433  << std::setw(8) << layerId << "\t" << std::setw(12) << scId << "\t" << std::setw(15)
434  << simClusters[scId].energy() << "\t" << std::setw(15) << SCenergy << "\t" << std::setw(14) << SCNumberOfHits
435  << "\t" << std::setw(18) << lcWithMaxEnergyInSC << "\t" << std::setw(15) << maxEnergyLCinSC << "\t"
436  << std::setw(20) << SCEnergyFractionInLC << "\n";
437 #endif
438  // Compute the correct normalization. Observe that this is the sum of the squares.
439  float invSCEnergyWeight = 0.f;
440  for (auto const& haf : lcsInSimCluster[scId][layerId].hits_and_fractions) {
441  invSCEnergyWeight += std::pow(haf.second * hitMap_->at(haf.first)->energy(), 2);
442  }
443  invSCEnergyWeight = 1.f / invSCEnergyWeight;
444  for (unsigned int i = 0; i < SCNumberOfHits; ++i) {
445  auto& sc_hitDetId = lcsInSimCluster[scId][layerId].hits_and_fractions[i].first;
446  auto& scFraction = lcsInSimCluster[scId][layerId].hits_and_fractions[i].second;
447 
448  bool hitWithLC = false;
449  if (scFraction == 0.f)
450  continue; //hopefully this should never happen
451  auto hit_find_in_LC = detIdToLayerClusterId_Map.find(sc_hitDetId);
452  if (hit_find_in_LC != detIdToLayerClusterId_Map.end())
453  hitWithLC = true;
454  auto itcheck = hitMap_->find(sc_hitDetId);
455  const HGCRecHit* hit = itcheck->second;
456  float hitEnergyWeight = hit->energy() * hit->energy();
457  for (auto& lcPair : lcsInSimCluster[scId][layerId].layerClusterIdToEnergyAndScore) {
458  unsigned int layerClusterId = lcPair.first;
459  float lcFraction = 0.f;
460 
461  if (hitWithLC) {
462  auto findHitIt = std::find(detIdToLayerClusterId_Map[sc_hitDetId].begin(),
463  detIdToLayerClusterId_Map[sc_hitDetId].end(),
464  hgcal::detIdInfoInCluster{layerClusterId, 0.f});
465  if (findHitIt != detIdToLayerClusterId_Map[sc_hitDetId].end())
466  lcFraction = findHitIt->fraction;
467  }
468  lcPair.second.second += std::min(std::pow(lcFraction - scFraction, 2), std::pow(scFraction, 2)) *
469  hitEnergyWeight * invSCEnergyWeight;
470 #ifdef EDM_ML_DEBUG
471  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
472  << "scDetId:\t" << (uint32_t)sc_hitDetId << "\tlayerClusterId:\t" << layerClusterId << "\t"
473  << "lcfraction,scfraction:\t" << lcFraction << ", " << scFraction << "\t"
474  << "hitEnergyWeight:\t" << hitEnergyWeight << "\t"
475  << "current score:\t" << lcPair.second.second << "\t"
476  << "invSCEnergyWeight:\t" << invSCEnergyWeight << "\n";
477 #endif
478  } // End of loop over LayerClusters linked to hits of this SimCluster
479  } // End of loop over hits of SimCluster on a Layer
480 #ifdef EDM_ML_DEBUG
481  if (lcsInSimCluster[scId][layerId].layerClusterIdToEnergyAndScore.empty())
482  LogDebug("LCToSCAssociatorByEnergyScoreImpl") << "SC Id: \t" << scId << "\tLC id:\t-1 "
483  << "\t score \t-1"
484  << "\n";
485 
486  for (const auto& lcPair : lcsInSimCluster[scId][layerId].layerClusterIdToEnergyAndScore) {
487  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
488  << "SC Id: \t" << scId << "\t LC id: \t" << lcPair.first << "\t score \t" << lcPair.second.second
489  << "\t shared energy:\t" << lcPair.second.first << "\t shared energy fraction:\t"
490  << (lcPair.second.first / SCenergy) << "\n";
491  }
492 #endif
493  } // End of loop over layers
494  } // End of loop over SimClusters
495 
496  return {scsInLayerCluster, lcsInSimCluster};
497 }
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
T const * product() const
Definition: Handle.h:70
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< std::vector< hgcal::simClusterOnLayer > > simClusterToLayerCluster
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
def unique(seq, keepstr=True)
Definition: tier0.py:24
double f[11][100]
const std::unordered_map< DetId, const HGCRecHit * > * hitMap_
Definition: DetId.h:17
std::shared_ptr< hgcal::RecHitTools > recHitTools_
std::vector< std::vector< std::pair< unsigned int, float > > > layerClusterToSimCluster
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
#define LogDebug(id)

Member Data Documentation

◆ hardScatterOnly_

const bool LCToSCAssociatorByEnergyScoreImpl::hardScatterOnly_
private

Definition at line 78 of file LCToSCAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ hitMap_

const std::unordered_map<DetId, const HGCRecHit *>* LCToSCAssociatorByEnergyScoreImpl::hitMap_
private

Definition at line 80 of file LCToSCAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ layers_

unsigned LCToSCAssociatorByEnergyScoreImpl::layers_
private

◆ productGetter_

edm::EDProductGetter const* LCToSCAssociatorByEnergyScoreImpl::productGetter_
private

Definition at line 82 of file LCToSCAssociatorByEnergyScoreImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

◆ recHitTools_

std::shared_ptr<hgcal::RecHitTools> LCToSCAssociatorByEnergyScoreImpl::recHitTools_
private