CMS 3D CMS Logo

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

#include <TSToSCAssociatorByEnergyScoreImpl.h>

Inheritance diagram for TSToSCAssociatorByEnergyScoreImpl:
hgcal::TracksterToSimClusterAssociatorBaseImpl

Public Member Functions

hgcal::RecoToSimCollectionTracksters associateRecoToSim (const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, const edm::Handle< SimClusterCollection > &sCCH) const override
 Associate a Trackster to SimClusters. More...
 
hgcal::SimToRecoCollectionTracksters associateSimToReco (const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, const edm::Handle< SimClusterCollection > &sCCH) const override
 Associate a SimCluster to Tracksters. More...
 
 TSToSCAssociatorByEnergyScoreImpl (edm::EDProductGetter const &, bool, std::shared_ptr< hgcal::RecHitTools >, const std::unordered_map< DetId, const HGCRecHit *> *)
 
- Public Member Functions inherited from hgcal::TracksterToSimClusterAssociatorBaseImpl
 TracksterToSimClusterAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~TracksterToSimClusterAssociatorBaseImpl ()
 Destructor. More...
 

Private Member Functions

hgcal::association makeConnections (const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, 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 40 of file TSToSCAssociatorByEnergyScoreImpl.h.

Constructor & Destructor Documentation

◆ TSToSCAssociatorByEnergyScoreImpl()

TSToSCAssociatorByEnergyScoreImpl::TSToSCAssociatorByEnergyScoreImpl ( 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 TSToSCAssociatorByEnergyScoreImpl.cc.

References layers_, and recHitTools_.

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

Member Function Documentation

◆ associateRecoToSim()

hgcal::RecoToSimCollectionTracksters TSToSCAssociatorByEnergyScoreImpl::associateRecoToSim ( const edm::Handle< ticl::TracksterCollection > &  tCH,
const edm::Handle< reco::CaloClusterCollection > &  lCCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
overridevirtual

Associate a Trackster to SimClusters.

Reimplemented from hgcal::TracksterToSimClusterAssociatorBaseImpl.

Definition at line 483 of file TSToSCAssociatorByEnergyScoreImpl.cc.

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

486  {
488  const auto& links = makeConnections(tCH, lCCH, sCCH);
489 
490  const auto& scsInTrackster = std::get<0>(links);
491  for (size_t tsId = 0; tsId < scsInTrackster.size(); ++tsId) {
492  for (auto& scPair : scsInTrackster[tsId]) {
493  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
494  << "Trackster Id:\t" << tsId << "\tSimCluster id:\t" << scPair.first << "\tscore:\t" << scPair.second << "\n";
495  // Fill AssociationMap
496  returnValue.insert(edm::Ref<ticl::TracksterCollection>(tCH, tsId), // Ref to TS
497  std::make_pair(edm::Ref<SimClusterCollection>(sCCH, scPair.first),
498  scPair.second) // Pair <Ref to SC, score>
499  );
500  }
501  }
502  return returnValue;
503 }
hgcal::association makeConnections(const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, const edm::Handle< SimClusterCollection > &sCCH) const
#define LogDebug(id)

◆ associateSimToReco()

hgcal::SimToRecoCollectionTracksters TSToSCAssociatorByEnergyScoreImpl::associateSimToReco ( const edm::Handle< ticl::TracksterCollection > &  tCH,
const edm::Handle< reco::CaloClusterCollection > &  lCCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
overridevirtual

Associate a SimCluster to Tracksters.

Reimplemented from hgcal::TracksterToSimClusterAssociatorBaseImpl.

Definition at line 505 of file TSToSCAssociatorByEnergyScoreImpl.cc.

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

508  {
510  const auto& links = makeConnections(tCH, lCCH, sCCH);
511  const auto& tssInSimCluster = std::get<1>(links);
512  for (size_t scId = 0; scId < tssInSimCluster.size(); ++scId) {
513  for (auto& tsPair : tssInSimCluster[scId].tracksterIdToEnergyAndScore) {
514  returnValue.insert(
515  edm::Ref<SimClusterCollection>(sCCH, scId), // Ref to SC
516  std::make_pair(edm::Ref<ticl::TracksterCollection>(tCH, tsPair.first), // Pair <Ref to TS,
517  std::make_pair(tsPair.second.first, tsPair.second.second)) // pair <energy, score> >
518  );
519  }
520  }
521  return returnValue;
522 }
hgcal::association makeConnections(const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, const edm::Handle< SimClusterCollection > &sCCH) const

◆ makeConnections()

hgcal::association TSToSCAssociatorByEnergyScoreImpl::makeConnections ( const edm::Handle< ticl::TracksterCollection > &  tCH,
const edm::Handle< reco::CaloClusterCollection > &  lCCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
private

Definition at line 15 of file TSToSCAssociatorByEnergyScoreImpl.cc.

References c, relativeConstraints::empty, mps_fire::end, edmPickEvents::event, f, spr::find(), cms::cuda::for(), h, hardScatterOnly_, hitMap_, mps_fire::i, dqmdumpme::last, hltEgammaHGCALIDVarsL1Seeded_cfi::layerClusters, LogDebug, or, funct::pow(), edm::Handle< T >::product(), jetsAK4_CHS_cff::sort, tier0::unique(), and AlignmentTracksFromVertexSelector_cfi::vertices.

Referenced by associateRecoToSim(), and associateSimToReco().

18  {
19  // Get collections
20  const auto& tracksters = *tCH.product();
21  const auto& layerClusters = *lCCH.product();
22  const auto& simClusters = *sCCH.product();
23  auto nTracksters = tracksters.size();
24 
25  //There shouldn't be any SimTracks from different crossings, but maybe they will be added later.
26  //At the moment there should be one SimTrack in each SimCluster.
27  auto nSimClusters = simClusters.size();
28  std::vector<size_t> sCIndices;
29  for (unsigned int scId = 0; scId < nSimClusters; ++scId) {
30  if (hardScatterOnly_ && (simClusters[scId].g4Tracks()[0].eventId().event() != 0 or
31  simClusters[scId].g4Tracks()[0].eventId().bunchCrossing() != 0)) {
32  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
33  << "Excluding SimCluster from event: " << simClusters[scId].g4Tracks()[0].eventId().event()
34  << " with BX: " << simClusters[scId].g4Tracks()[0].eventId().bunchCrossing() << std::endl;
35  continue;
36  }
37  sCIndices.emplace_back(scId);
38  }
39  nSimClusters = sCIndices.size();
40 
41  // Initialize tssInSimCluster. To be returned outside, since it contains the
42  // information to compute the SimCluster-To-Trackster score.
43  // tssInSimCluster[scId]:
44  hgcal::simClusterToTrackster tssInSimCluster;
45  tssInSimCluster.resize(nSimClusters);
46  for (unsigned int i = 0; i < nSimClusters; ++i) {
47  tssInSimCluster[i].simClusterId = i;
48  tssInSimCluster[i].energy = 0.f;
49  tssInSimCluster[i].hits_and_fractions.clear();
50  }
51 
52  // Fill detIdToSimClusterId_Map and update tssInSimCluster
53  std::unordered_map<DetId, std::vector<hgcal::detIdInfoInCluster>> detIdToSimClusterId_Map;
54  for (const auto& scId : sCIndices) {
55  const auto& hits_and_fractions = simClusters[scId].hits_and_fractions();
56  for (const auto& it_haf : hits_and_fractions) {
57  const auto hitid = it_haf.first;
58  const auto itcheck = hitMap_->find(hitid);
59  if (itcheck != hitMap_->end()) {
60  const auto hit_find_it = detIdToSimClusterId_Map.find(hitid);
61  if (hit_find_it == detIdToSimClusterId_Map.end()) {
62  detIdToSimClusterId_Map[hitid] = std::vector<hgcal::detIdInfoInCluster>();
63  }
64  detIdToSimClusterId_Map[hitid].emplace_back(scId, it_haf.second);
65 
66  const HGCRecHit* hit = itcheck->second;
67  tssInSimCluster[scId].energy += it_haf.second * hit->energy();
68  tssInSimCluster[scId].hits_and_fractions.emplace_back(hitid, it_haf.second);
69  }
70  }
71  } // end of loop over SimClusters
72 
73 #ifdef EDM_ML_DEBUG
74  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
75  << "tssInSimCluster INFO (Only SimCluster filled at the moment)" << std::endl;
76  for (size_t sc = 0; sc < tssInSimCluster.size(); ++sc) {
77  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "For SimCluster Idx: " << sc << " we have: " << std::endl;
78  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
79  << "\tSimClusterIdx:\t" << tssInSimCluster[sc].simClusterId << std::endl;
80  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "\tEnergy:\t" << tssInSimCluster[sc].energy << std::endl;
81  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "\t# of clusters:\t" << layerClusters.size() << std::endl;
82  double tot_energy = 0.;
83  for (auto const& haf : tssInSimCluster[sc].hits_and_fractions) {
84  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
85  << "\tHits/fraction/energy: " << (uint32_t)haf.first << "/" << haf.second << "/"
86  << haf.second * hitMap_->at(haf.first)->energy() << std::endl;
87  tot_energy += haf.second * hitMap_->at(haf.first)->energy();
88  }
89  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "\tTot Sum haf: " << tot_energy << std::endl;
90  for (auto const& ts : tssInSimCluster[sc].tracksterIdToEnergyAndScore) {
91  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
92  << "\ttsIdx/energy/score: " << ts.first << "/" << ts.second.first << "/" << ts.second.second << std::endl;
93  }
94  }
95 
96  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "detIdToSimClusterId_Map INFO" << std::endl;
97  for (auto const& detId : detIdToSimClusterId_Map) {
98  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
99  << "For detId: " << (uint32_t)detId.first
100  << " we have found the following connections with SimClusters:" << std::endl;
101  for (auto const& sc : detId.second) {
102  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
103  << "\tSimCluster Id: " << sc.clusterId << " with fraction: " << sc.fraction
104  << " and energy: " << sc.fraction * hitMap_->at(detId.first)->energy() << std::endl;
105  }
106  }
107 #endif
108 
109  // Fill detIdToLayerClusterId_Map and scsInTrackster; update tssInSimCluster
110  std::unordered_map<DetId, std::vector<hgcal::detIdInfoInCluster>> detIdToLayerClusterId_Map;
111  // this contains the ids of the simclusters contributing with at least one
112  // hit to the Trackster. To be returned since this contains the information
113  // to compute the Trackster-To-SimCluster score.
114  hgcal::tracksterToSimCluster scsInTrackster; //[tsId][scId]->(energy,score)
115  scsInTrackster.resize(nTracksters);
116 
117  for (unsigned int tsId = 0; tsId < nTracksters; ++tsId) {
118  for (unsigned int i = 0; i < tracksters[tsId].vertices().size(); ++i) {
119  const auto lcId = tracksters[tsId].vertices(i);
120  const auto lcFractionInTs = 1.f / tracksters[tsId].vertex_multiplicity(i);
121 
122  const std::vector<std::pair<DetId, float>>& hits_and_fractions = layerClusters[lcId].hitsAndFractions();
123  unsigned int numberOfHitsInLC = hits_and_fractions.size();
124 
125  for (unsigned int hitId = 0; hitId < numberOfHitsInLC; hitId++) {
126  const auto rh_detid = hits_and_fractions[hitId].first;
127  const auto rhFraction = hits_and_fractions[hitId].second;
128 
129  const auto hit_find_in_LC = detIdToLayerClusterId_Map.find(rh_detid);
130  if (hit_find_in_LC == detIdToLayerClusterId_Map.end()) {
131  detIdToLayerClusterId_Map[rh_detid] = std::vector<hgcal::detIdInfoInCluster>();
132  }
133  detIdToLayerClusterId_Map[rh_detid].emplace_back(lcId, rhFraction);
134 
135  const auto hit_find_in_SC = detIdToSimClusterId_Map.find(rh_detid);
136 
137  if (hit_find_in_SC != detIdToSimClusterId_Map.end()) {
138  const auto itcheck = hitMap_->find(rh_detid);
139  const HGCRecHit* hit = itcheck->second;
140  //Loops through all the simclusters that have the layer cluster rechit under study
141  //Here is time to update the tssInSimCluster and connect the SimCluster with all
142  //the Tracksters that have the current rechit detid matched.
143  for (const auto& h : hit_find_in_SC->second) {
144  //tssInSimCluster[simclusterId][layerclusterId]-> (energy,score)
145  //SC_i - > TS_j, TS_k, ...
146  tssInSimCluster[h.clusterId].tracksterIdToEnergyAndScore[tsId].first +=
147  lcFractionInTs * h.fraction * hit->energy();
148  //TS_i -> SC_j, SC_k, ...
149  scsInTrackster[tsId].emplace_back(h.clusterId, 0.f);
150  }
151  }
152  } // End loop over hits on a LayerCluster
153  } // End loop over LayerClusters in Trackster
154  } // End of loop over Tracksters
155 
156 #ifdef EDM_ML_DEBUG
157  for (unsigned int tsId = 0; tsId < nTracksters; ++tsId) {
158  for (const auto& lcId : tracksters[tsId].vertices()) {
159  const auto& hits_and_fractions = layerClusters[lcId].hitsAndFractions();
160  unsigned int numberOfHitsInLC = hits_and_fractions.size();
161 
162  // This vector will store, for each hit in the Layercluster, the index of
163  // the SimCluster that contributed the most, in terms of energy, to it.
164  // Special values are:
165  //
166  // -2 --> the reconstruction fraction of the RecHit is 0 (used in the past to monitor Halo Hits)
167  // -3 --> same as before with the added condition that no SimCluster has been linked to this RecHit
168  // -1 --> the reco fraction is >0, but no SimCluster has been linked to it
169  // >=0 --> index of the linked SimCluster
170  std::vector<int> hitsToSimClusterId(numberOfHitsInLC);
171  // This will store the index of the SimCluster linked to the LayerCluster that has the largest number of hits in common.
172  int maxSCId_byNumberOfHits = -1;
173  // This will store the maximum number of shared hits between a LayerCluster and a SimCluster
174  unsigned int maxSCNumberOfHitsInLC = 0;
175  // This will store the index of the SimCluster linked to the LayerCluster that has the largest energy in common.
176  int maxSCId_byEnergy = -1;
177  // This will store the maximum number of shared energy between a LayerCluster and a SimCluster
178  float maxEnergySharedLCandSC = 0.f;
179  // This will store the fraction of the LayerCluster energy shared with the best(energy) SimCluster: e_shared/lc_energy
180  float energyFractionOfLCinSC = 0.f;
181  // This will store the fraction of the SimCluster energy shared with the Trackster: e_shared/sc_energy
182  float energyFractionOfSCinLC = 0.f;
183  std::unordered_map<unsigned, unsigned> occurrencesSCinLC;
184  unsigned int numberOfNoiseHitsInLC = 0;
185  std::unordered_map<unsigned, float> SCEnergyInLC;
186 
187  for (unsigned int hitId = 0; hitId < numberOfHitsInLC; hitId++) {
188  const auto rh_detid = hits_and_fractions[hitId].first;
189  const auto rhFraction = hits_and_fractions[hitId].second;
190 
191  const auto hit_find_in_SC = detIdToSimClusterId_Map.find(rh_detid);
192 
193  // if the fraction is zero or the hit does not belong to any SimCluster,
194  // set the SimCluster Id for the hit to -1; this will
195  // contribute to the number of noise hits
196 
197  // MR Remove the case in which the fraction is 0, since this could be a
198  // real hit that has been marked as halo.
199  if (rhFraction == 0.) {
200  hitsToSimClusterId[hitId] = -2;
201  }
202  //Now check if there are SimClusters linked to this rechit of the layercluster
203  if (hit_find_in_SC == detIdToSimClusterId_Map.end()) {
204  hitsToSimClusterId[hitId] -= 1;
205  } else {
206  const auto itcheck = hitMap_->find(rh_detid);
207  const HGCRecHit* hit = itcheck->second;
208  auto maxSCEnergyInLC = 0.f;
209  auto maxSCId = -1;
210  //Loop through all the linked SimClusters
211  for (const auto& h : hit_find_in_SC->second) {
212  SCEnergyInLC[h.clusterId] += h.fraction * hit->energy();
213  // Keep track of which SimCluster contributed the most, in terms
214  // of energy, to this specific Layer Cluster.
215  if (SCEnergyInLC[h.clusterId] > maxSCEnergyInLC) {
216  maxSCEnergyInLC = SCEnergyInLC[h.clusterId];
217  maxSCId = h.clusterId;
218  }
219  }
220  hitsToSimClusterId[hitId] = maxSCId;
221  }
222  } // End loop over hits on a LayerCluster
223 
224  for (const auto& c : hitsToSimClusterId) {
225  if (c < 0) {
226  numberOfNoiseHitsInLC++;
227  } else {
228  occurrencesSCinLC[c]++;
229  }
230  }
231 
232  for (const auto& c : occurrencesSCinLC) {
233  if (c.second > maxSCNumberOfHitsInLC) {
234  maxSCId_byNumberOfHits = c.first;
235  maxSCNumberOfHitsInLC = c.second;
236  }
237  }
238 
239  for (const auto& c : SCEnergyInLC) {
240  if (c.second > maxEnergySharedLCandSC) {
241  maxSCId_byEnergy = c.first;
242  maxEnergySharedLCandSC = c.second;
243  }
244  }
245 
246  float totalSCEnergyOnLayer = 0.f;
247  if (maxSCId_byEnergy >= 0) {
248  totalSCEnergyOnLayer = tssInSimCluster[maxSCId_byEnergy].energy;
249  energyFractionOfSCinLC = maxEnergySharedLCandSC / totalSCEnergyOnLayer;
250  if (tracksters[tsId].raw_energy() > 0.f) {
251  energyFractionOfLCinSC = maxEnergySharedLCandSC / tracksters[tsId].raw_energy();
252  }
253  }
254 
255  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
256  << std::setw(12) << "TracksterID:\t" << std::setw(12) << "layerCluster\t" << std::setw(10) << "lc energy\t"
257  << std::setw(5) << "nhits\t" << std::setw(12) << "noise hits\t" << std::setw(22) << "maxSCId_byNumberOfHits\t"
258  << std::setw(8) << "nhitsSC\t" << std::setw(13) << "maxSCId_byEnergy\t" << std::setw(20)
259  << "maxEnergySharedLCandSC\t" << std::setw(22) << "totalSCEnergyOnLayer\t" << std::setw(22)
260  << "energyFractionOfLCinSC\t" << std::setw(25) << "energyFractionOfSCinLC\t"
261  << "\n";
262  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
263  << std::setw(12) << tsId << "\t" << std::setw(12) << lcId << "\t" << std::setw(10)
264  << tracksters[tsId].raw_energy() << "\t" << std::setw(5) << numberOfHitsInLC << "\t" << std::setw(12)
265  << numberOfNoiseHitsInLC << "\t" << std::setw(22) << maxSCId_byNumberOfHits << "\t" << std::setw(8)
266  << maxSCNumberOfHitsInLC << "\t" << std::setw(13) << maxSCId_byEnergy << "\t" << std::setw(20)
267  << maxEnergySharedLCandSC << "\t" << std::setw(22) << totalSCEnergyOnLayer << "\t" << std::setw(22)
268  << energyFractionOfLCinSC << "\t" << std::setw(25) << energyFractionOfSCinLC << "\n";
269  } // End of loop over LayerClusters in Trackster
270  } // End of loop over Tracksters
271 
272  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
273  << "Improved tssInSimCluster INFO (Now containing the linked tracksters id and energy - score still empty)"
274  << std::endl;
275  for (size_t sc = 0; sc < tssInSimCluster.size(); ++sc) {
276  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "For SimCluster Idx: " << sc << " we have: " << std::endl;
277  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
278  << " SimClusterIdx: " << tssInSimCluster[sc].simClusterId << std::endl;
279  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "\tEnergy:\t" << tssInSimCluster[sc].energy << std::endl;
280  double tot_energy = 0.;
281  for (auto const& haf : tssInSimCluster[sc].hits_and_fractions) {
282  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
283  << "\tHits/fraction/energy: " << (uint32_t)haf.first << "/" << haf.second << "/"
284  << haf.second * hitMap_->at(haf.first)->energy() << std::endl;
285  tot_energy += haf.second * hitMap_->at(haf.first)->energy();
286  }
287  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "\tTot Sum haf: " << tot_energy << std::endl;
288  for (auto const& ts : tssInSimCluster[sc].tracksterIdToEnergyAndScore) {
289  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
290  << "\ttsIdx/energy/score: " << ts.first << "/" << ts.second.first << "/" << ts.second.second << std::endl;
291  }
292  }
293 
294  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "Improved detIdToSimClusterId_Map INFO" << std::endl;
295  for (auto const& sc : detIdToSimClusterId_Map) {
296  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
297  << "For detId: " << (uint32_t)sc.first
298  << " we have found the following connections with SimClusters:" << std::endl;
299  for (auto const& sclu : sc.second) {
300  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
301  << " SimCluster Id: " << sclu.clusterId << " with fraction: " << sclu.fraction
302  << " and energy: " << sclu.fraction * hitMap_->at(sc.first)->energy() << std::endl;
303  }
304  }
305 #endif
306 
307  // Update scsInTrackster; compute the score Trackster-to-SimCluster,
308  // together with the returned AssociationMap
309  for (unsigned int tsId = 0; tsId < nTracksters; ++tsId) {
310  // The SimClusters contributing to the Trackster's LayerClusters should already be unique.
311  // find the unique SimClusters id contributing to the Trackster's LayerClusters
312  std::sort(scsInTrackster[tsId].begin(), scsInTrackster[tsId].end());
313  auto last = std::unique(scsInTrackster[tsId].begin(), scsInTrackster[tsId].end());
314  scsInTrackster[tsId].erase(last, scsInTrackster[tsId].end());
315 
316  // If a reconstructed Trackster has energy 0 but is linked to a
317  // SimCluster, assigned score 1
318  if (tracksters[tsId].raw_energy() == 0. && !scsInTrackster[tsId].empty()) {
319  for (auto& scPair : scsInTrackster[tsId]) {
320  scPair.second = 1.;
321  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
322  << "TracksterId:\t " << tsId << "\tSC id:\t" << scPair.first << "\tscore\t " << scPair.second << "\n";
323  }
324  continue;
325  }
326 
327  float invTracksterEnergyWeight = 0.f;
328  for (unsigned int i = 0; i < tracksters[tsId].vertices().size(); ++i) {
329  const auto lcId = tracksters[tsId].vertices(i);
330  const auto lcFractionInTs = 1.f / tracksters[tsId].vertex_multiplicity(i);
331 
332  const auto& hits_and_fractions = layerClusters[lcId].hitsAndFractions();
333  // Compute the correct normalization
334  for (auto const& haf : hits_and_fractions) {
335  invTracksterEnergyWeight += (lcFractionInTs * haf.second * hitMap_->at(haf.first)->energy()) *
336  (lcFractionInTs * haf.second * hitMap_->at(haf.first)->energy());
337  }
338  }
339  invTracksterEnergyWeight = 1.f / invTracksterEnergyWeight;
340 
341  for (unsigned int i = 0; i < tracksters[tsId].vertices().size(); ++i) {
342  const auto lcId = tracksters[tsId].vertices(i);
343  const auto lcFractionInTs = 1.f / tracksters[tsId].vertex_multiplicity(i);
344 
345  const auto& hits_and_fractions = layerClusters[lcId].hitsAndFractions();
346  unsigned int numberOfHitsInLC = hits_and_fractions.size();
347  for (unsigned int i = 0; i < numberOfHitsInLC; ++i) {
348  DetId rh_detid = hits_and_fractions[i].first;
349  float rhFraction = hits_and_fractions[i].second * lcFractionInTs;
350 
351  const bool hitWithSC = (detIdToSimClusterId_Map.find(rh_detid) != detIdToSimClusterId_Map.end());
352 
353  const auto itcheck = hitMap_->find(rh_detid);
354  const HGCRecHit* hit = itcheck->second;
355  float hitEnergyWeight = hit->energy() * hit->energy();
356 
357  for (auto& scPair : scsInTrackster[tsId]) {
358  float scFraction = 0.f;
359  if (hitWithSC) {
360  const auto findHitIt = std::find(detIdToSimClusterId_Map[rh_detid].begin(),
361  detIdToSimClusterId_Map[rh_detid].end(),
362  hgcal::detIdInfoInCluster{scPair.first, 0.f});
363  if (findHitIt != detIdToSimClusterId_Map[rh_detid].end())
364  scFraction = findHitIt->fraction;
365  }
366  scPair.second +=
367  (rhFraction - scFraction) * (rhFraction - scFraction) * hitEnergyWeight * invTracksterEnergyWeight;
368 #ifdef EDM_ML_DEBUG
369  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
370  << "rh_detid:\t" << (uint32_t)rh_detid << "\ttracksterId:\t" << tsId << "\t"
371  << "rhfraction,scFraction:\t" << rhFraction << ", " << scFraction << "\t"
372  << "hitEnergyWeight:\t" << hitEnergyWeight << "\t"
373  << "current score:\t" << scPair.second << "\t"
374  << "invTracksterEnergyWeight:\t" << invTracksterEnergyWeight << "\n";
375 #endif
376  }
377  } // End of loop over Hits within a LayerCluster
378  } // End of loop over LayerClusters in Trackster
379 
380 #ifdef EDM_ML_DEBUG
381  if (scsInTrackster[tsId].empty())
382  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "trackster Id:\t" << tsId << "\tSC id:\t-1"
383  << "\tscore\t-1\n";
384 #endif
385  } // End of loop over Tracksters
386 
387  // Compute the SimCluster-To-Trackster score
388  for (const auto& scId : sCIndices) {
389  float invSCEnergyWeight = 0.f;
390 
391  const unsigned int SCNumberOfHits = tssInSimCluster[scId].hits_and_fractions.size();
392  if (SCNumberOfHits == 0)
393  continue;
394 #ifdef EDM_ML_DEBUG
395  int tsWithMaxEnergyInSC = -1;
396  //energy of the most energetic TS from all that were linked to SC
397  float maxEnergyTSinSC = 0.f;
398  float SCenergy = tssInSimCluster[scId].energy;
399  //most energetic TS from all TSs linked to SC over SC energy.
400  float SCEnergyFractionInTS = 0.f;
401  for (const auto& ts : tssInSimCluster[scId].tracksterIdToEnergyAndScore) {
402  if (ts.second.first > maxEnergyTSinSC) {
403  maxEnergyTSinSC = ts.second.first;
404  tsWithMaxEnergyInSC = ts.first;
405  }
406  }
407  if (SCenergy > 0.f)
408  SCEnergyFractionInTS = maxEnergyTSinSC / SCenergy;
409 
410  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
411  << std::setw(12) << "simcluster\t" << std::setw(15) << "sc total energy\t" << std::setw(15)
412  << "scEnergyOnLayer\t" << std::setw(14) << "SCNhitsOnLayer\t" << std::setw(18) << "tsWithMaxEnergyInSC\t"
413  << std::setw(15) << "maxEnergyTSinSC\t" << std::setw(20) << "SCEnergyFractionInTS"
414  << "\n";
415  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
416  << std::setw(12) << scId << "\t" << std::setw(15) << simClusters[scId].energy() << "\t" << std::setw(15)
417  << SCenergy << "\t" << std::setw(14) << SCNumberOfHits << "\t" << std::setw(18) << tsWithMaxEnergyInSC << "\t"
418  << std::setw(15) << maxEnergyTSinSC << "\t" << std::setw(20) << SCEnergyFractionInTS << "\n";
419 #endif
420  // Compute the correct normalization
421  for (auto const& haf : tssInSimCluster[scId].hits_and_fractions) {
422  invSCEnergyWeight += std::pow(haf.second * hitMap_->at(haf.first)->energy(), 2);
423  }
424  invSCEnergyWeight = 1.f / invSCEnergyWeight;
425 
426  for (unsigned int i = 0; i < SCNumberOfHits; ++i) {
427  auto& sc_hitDetId = tssInSimCluster[scId].hits_and_fractions[i].first;
428  auto& scFraction = tssInSimCluster[scId].hits_and_fractions[i].second;
429 
430  bool hitWithLC = false;
431  if (scFraction == 0.f)
432  continue; // hopefully this should never happen
433  const auto hit_find_in_LC = detIdToLayerClusterId_Map.find(sc_hitDetId);
434  if (hit_find_in_LC != detIdToLayerClusterId_Map.end())
435  hitWithLC = true;
436  const auto itcheck = hitMap_->find(sc_hitDetId);
437  const HGCRecHit* hit = itcheck->second;
438  float hitEnergyWeight = hit->energy() * hit->energy();
439  for (auto& tsPair : tssInSimCluster[scId].tracksterIdToEnergyAndScore) {
440  unsigned int tsId = tsPair.first;
441  float tsFraction = 0.f;
442 
443  for (unsigned int i = 0; i < tracksters[tsId].vertices().size(); ++i) {
444  const auto lcId = tracksters[tsId].vertices(i);
445  const auto lcFractionInTs = 1.f / tracksters[tsId].vertex_multiplicity(i);
446 
447  if (hitWithLC) {
448  const auto findHitIt = std::find(detIdToLayerClusterId_Map[sc_hitDetId].begin(),
449  detIdToLayerClusterId_Map[sc_hitDetId].end(),
450  hgcal::detIdInfoInCluster{lcId, 0.f});
451  if (findHitIt != detIdToLayerClusterId_Map[sc_hitDetId].end())
452  tsFraction = findHitIt->fraction * lcFractionInTs;
453  }
454  tsPair.second.second +=
455  (tsFraction - scFraction) * (tsFraction - scFraction) * hitEnergyWeight * invSCEnergyWeight;
456 #ifdef EDM_ML_DEBUG
457  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
458  << "SCDetId:\t" << (uint32_t)sc_hitDetId << "\tTracksterId:\t" << tsId << "\t"
459  << "tsFraction, scFraction:\t" << tsFraction << ", " << scFraction << "\t"
460  << "hitEnergyWeight:\t" << hitEnergyWeight << "\t"
461  << "current score:\t" << tsPair.second.second << "\t"
462  << "invSCEnergyWeight:\t" << invSCEnergyWeight << "\n";
463 #endif
464  } // End of loop over Trackster's LayerClusters
465  } // End of loop over Tracksters linked to hits of this SimCluster
466  } // End of loop over hits of SimCluster on a Layer
467 #ifdef EDM_ML_DEBUG
468  if (tssInSimCluster[scId].tracksterIdToEnergyAndScore.empty())
469  LogDebug("TSToSCAssociatorByEnergyScoreImpl") << "SC Id:\t" << scId << "\tTS id:\t-1 "
470  << "\tscore\t-1\n";
471 
472  for (const auto& tsPair : tssInSimCluster[scId].tracksterIdToEnergyAndScore) {
473  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
474  << "SC Id: \t" << scId << "\t TS id: \t" << tsPair.first << "\t score \t" << tsPair.second.second
475  << "\t shared energy:\t" << tsPair.second.first << "\t shared energy fraction:\t"
476  << (tsPair.second.first / SCenergy) << "\n";
477  }
478 #endif
479  } // End loop over SimCluster indices
480  return {scsInTrackster, tssInSimCluster};
481 }
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
const std::unordered_map< DetId, const HGCRecHit * > * hitMap_
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
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]
std::vector< std::vector< std::pair< unsigned int, float > > > tracksterToSimCluster
std::vector< hgcal::simClusterOnBLayer > simClusterToTrackster
Definition: DetId.h:17
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 TSToSCAssociatorByEnergyScoreImpl::hardScatterOnly_
private

Definition at line 56 of file TSToSCAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ hitMap_

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

Definition at line 58 of file TSToSCAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ layers_

unsigned TSToSCAssociatorByEnergyScoreImpl::layers_
private

◆ productGetter_

edm::EDProductGetter const* TSToSCAssociatorByEnergyScoreImpl::productGetter_
private

Definition at line 60 of file TSToSCAssociatorByEnergyScoreImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

◆ recHitTools_

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