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:
ticl::TracksterToSimClusterAssociatorBaseImpl

Public Member Functions

ticl::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...
 
ticl::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 unsigned int > *, std::vector< const HGCRecHit *> &hits)
 
- Public Member Functions inherited from ticl::TracksterToSimClusterAssociatorBaseImpl
 TracksterToSimClusterAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~TracksterToSimClusterAssociatorBaseImpl ()
 Destructor. More...
 

Private Member Functions

ticl::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 unsigned int > * hitMap_
 
std::vector< const HGCRecHit * > hits_
 
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 unsigned int > *  hitMap,
std::vector< const HGCRecHit *> &  hits 
)
explicit

Definition at line 6 of file TSToSCAssociatorByEnergyScoreImpl.cc.

References layers_, and recHitTools_.

12  : hardScatterOnly_(hardScatterOnly),
13  recHitTools_(recHitTools),
14  hitMap_(hitMap),
15  hits_(hits),
17  layers_ = recHitTools_->lastLayerBH();
18 }
std::shared_ptr< hgcal::RecHitTools > recHitTools_
const std::unordered_map< DetId, const unsigned int > * hitMap_
EDProductGetter const * productGetter(std::atomic< void const *> const &iCache)

Member Function Documentation

◆ associateRecoToSim()

ticl::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 ticl::TracksterToSimClusterAssociatorBaseImpl.

Definition at line 492 of file TSToSCAssociatorByEnergyScoreImpl.cc.

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

495  {
497  const auto& links = makeConnections(tCH, lCCH, sCCH);
498 
499  const auto& scsInTrackster = std::get<0>(links);
500  for (size_t tsId = 0; tsId < scsInTrackster.size(); ++tsId) {
501  for (auto& scPair : scsInTrackster[tsId]) {
502  LogDebug("TSToSCAssociatorByEnergyScoreImpl")
503  << "Trackster Id:\t" << tsId << "\tSimCluster id:\t" << scPair.first << "\tscore:\t" << scPair.second << "\n";
504  // Fill AssociationMap
505  returnValue.insert(edm::Ref<ticl::TracksterCollection>(tCH, tsId), // Ref to TS
506  std::make_pair(edm::Ref<SimClusterCollection>(sCCH, scPair.first),
507  scPair.second) // Pair <Ref to SC, score>
508  );
509  }
510  }
511  return returnValue;
512 }
ticl::association makeConnections(const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, const edm::Handle< SimClusterCollection > &sCCH) const
#define LogDebug(id)

◆ associateSimToReco()

ticl::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 ticl::TracksterToSimClusterAssociatorBaseImpl.

Definition at line 514 of file TSToSCAssociatorByEnergyScoreImpl.cc.

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

517  {
519  const auto& links = makeConnections(tCH, lCCH, sCCH);
520  const auto& tssInSimCluster = std::get<1>(links);
521  for (size_t scId = 0; scId < tssInSimCluster.size(); ++scId) {
522  for (auto& tsPair : tssInSimCluster[scId].tracksterIdToEnergyAndScore) {
523  returnValue.insert(
524  edm::Ref<SimClusterCollection>(sCCH, scId), // Ref to SC
525  std::make_pair(edm::Ref<ticl::TracksterCollection>(tCH, tsPair.first), // Pair <Ref to TS,
526  std::make_pair(tsPair.second.first, tsPair.second.second)) // pair <energy, score> >
527  );
528  }
529  }
530  return returnValue;
531 }
ticl::association makeConnections(const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, const edm::Handle< SimClusterCollection > &sCCH) const

◆ makeConnections()

ticl::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 20 of file TSToSCAssociatorByEnergyScoreImpl.cc.

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

Referenced by associateRecoToSim(), and associateSimToReco().

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

Member Data Documentation

◆ hardScatterOnly_

const bool TSToSCAssociatorByEnergyScoreImpl::hardScatterOnly_
private

Definition at line 57 of file TSToSCAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ hitMap_

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

Definition at line 59 of file TSToSCAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ hits_

std::vector<const HGCRecHit *> TSToSCAssociatorByEnergyScoreImpl::hits_
private

Definition at line 60 of file TSToSCAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ layers_

unsigned TSToSCAssociatorByEnergyScoreImpl::layers_
private

◆ productGetter_

edm::EDProductGetter const* TSToSCAssociatorByEnergyScoreImpl::productGetter_
private

Definition at line 62 of file TSToSCAssociatorByEnergyScoreImpl.h.

Referenced by associateRecoToSim(), and associateSimToReco().

◆ recHitTools_

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