CMS 3D CMS Logo

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

#include <TSToSimTSAssociatorByEnergyScoreImpl.h>

Inheritance diagram for TSToSimTSAssociatorByEnergyScoreImpl:
hgcal::TracksterToSimTracksterAssociatorBaseImpl

Public Member Functions

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

Private Member Functions

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

Constructor & Destructor Documentation

◆ TSToSimTSAssociatorByEnergyScoreImpl()

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

Definition at line 5 of file TSToSimTSAssociatorByEnergyScoreImpl.cc.

References layers_, and recHitTools_.

10  : hardScatterOnly_(hardScatterOnly), recHitTools_(recHitTools), hitMap_(hitMap), productGetter_(&productGetter) {
11  layers_ = recHitTools_->lastLayerBH();
12 }
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::RecoToSimCollectionSimTracksters TSToSimTSAssociatorByEnergyScoreImpl::associateRecoToSim ( const edm::Handle< ticl::TracksterCollection > &  tCH,
const edm::Handle< reco::CaloClusterCollection > &  lCCH,
const edm::Handle< ticl::TracksterCollection > &  sTCH 
) const
overridevirtual

Associate a Trackster to SimClusters.

Reimplemented from hgcal::TracksterToSimTracksterAssociatorBaseImpl.

Definition at line 437 of file TSToSimTSAssociatorByEnergyScoreImpl.cc.

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

440  {
442  const auto& links = makeConnections(tCH, lCCH, sTCH);
443 
444  const auto& stsInTrackster = std::get<0>(links);
445  for (size_t tsId = 0; tsId < stsInTrackster.size(); ++tsId) {
446  for (auto& stPair : stsInTrackster[tsId]) {
447  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "Trackster Id:\t" << tsId << "\tSimTrackster id:\t"
448  << stPair.first << "\tscore:\t" << stPair.second.second << "\n";
449  // Fill AssociationMap
450  returnValue.insert(
451  edm::Ref<ticl::TracksterCollection>(tCH, tsId), // Ref to TS
452  std::make_pair(edm::Ref<ticl::TracksterCollection>(sTCH, stPair.first), //Pair <Refo to TS>
453  std::make_pair(stPair.second.first, stPair.second.second)) // Pair <energy, score>
454  );
455  }
456  }
457  return returnValue;
458 }
hgcal::association makeConnections(const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, const edm::Handle< ticl::TracksterCollection > &sTCH) const
#define LogDebug(id)

◆ associateSimToReco()

hgcal::SimToRecoCollectionSimTracksters TSToSimTSAssociatorByEnergyScoreImpl::associateSimToReco ( const edm::Handle< ticl::TracksterCollection > &  tCH,
const edm::Handle< reco::CaloClusterCollection > &  lCCH,
const edm::Handle< ticl::TracksterCollection > &  sTCH 
) const
overridevirtual

Associate a SimCluster to Tracksters.

Reimplemented from hgcal::TracksterToSimTracksterAssociatorBaseImpl.

Definition at line 460 of file TSToSimTSAssociatorByEnergyScoreImpl.cc.

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

463  {
465  const auto& links = makeConnections(tCH, lCCH, sTCH);
466  const auto& tssInSimTrackster = std::get<1>(links);
467  for (size_t stId = 0; stId < tssInSimTrackster.size(); ++stId) {
468  for (auto& tsPair : tssInSimTrackster[stId].tracksterIdToEnergyAndScore) {
469  returnValue.insert(
470  edm::Ref<ticl::TracksterCollection>(sTCH, stId), // Ref to ST
471  std::make_pair(edm::Ref<ticl::TracksterCollection>(tCH, tsPair.first), // Pair <Ref to TS,
472  std::make_pair(tsPair.second.first, tsPair.second.second)) // pair <energy, score> >
473  );
474  }
475  }
476  return returnValue;
477 }
hgcal::association makeConnections(const edm::Handle< ticl::TracksterCollection > &tCH, const edm::Handle< reco::CaloClusterCollection > &lCCH, const edm::Handle< ticl::TracksterCollection > &sTCH) const

◆ makeConnections()

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

Definition at line 14 of file TSToSimTSAssociatorByEnergyScoreImpl.cc.

References HltBtagPostValidation_cff::c, relativeConstraints::empty, mps_fire::end, hcalRecHitTable_cff::energy, f, spr::find(), cms::cuda::for(), HLT_2024v10_cff::fraction, hitMap_, mps_fire::i, dqmdumpme::last, hltEgammaHGCALIDVarsL1Seeded_cfi::layerClusters, LogDebug, funct::pow(), edm::Handle< T >::product(), jetUpdater_cfi::sort, tier0::unique(), and AlignmentTracksFromVertexSelector_cfi::vertices.

Referenced by associateRecoToSim(), and associateSimToReco().

17  {
18  // Get collections
19  const auto& tracksters = *tCH.product();
20  const auto& layerClusters = *lCCH.product();
21  const auto& simTracksters = *sTCH.product();
22  auto nTracksters = tracksters.size();
23 
24  //There shouldn't be any SimTrackster without vertices, but maybe they will be added later.
25  auto nSimTracksters = simTracksters.size();
26  std::vector<size_t> sTIndices;
27  for (unsigned int stId = 0; stId < nSimTracksters; ++stId) {
28  if (simTracksters[stId].vertices().empty()) {
29  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
30  << "Excluding SimTrackster " << stId << " witH no vertices!" << std::endl;
31  continue;
32  }
33  sTIndices.emplace_back(stId);
34  }
35  nSimTracksters = sTIndices.size();
36 
37  // Initialize tssInSimTrackster. It contains the simTracksterOnLayer structure for all CaloParticles in each layer and
38  // among other the information to compute the SimTrackster-To-Trackster score. It is one of the two objects that
39  // build the output of the makeConnections function.
40  // tssInSimTrackster[stId]
41  hgcal::simTracksterToTrackster tssInSimTrackster;
42  tssInSimTrackster.resize(nSimTracksters);
43  for (unsigned int i = 0; i < nSimTracksters; ++i) {
44  tssInSimTrackster[i].simTracksterId = i;
45  tssInSimTrackster[i].energy = 0.f;
46  tssInSimTrackster[i].lcs_and_fractions.clear();
47  }
48 
49  // Fill lcToSimTracksterId_Map and update tssInSimTrackster
50  // The lcToSimTracksterId_Map is used to connect a LayerCluster, via its id (key), with all the SimTracksters that
51  // contributed to that LayerCluster by storing the SimTrackster id and the fraction of the LayerCluster's energy
52  // in which the SimTrackster contributed.
53  std::unordered_map<int, std::vector<hgcal::lcInfoInTrackster>> lcToSimTracksterId_Map;
54  for (const auto& stId : sTIndices) {
55  const auto& lcs = simTracksters[stId].vertices();
56  int lcInSimTst = 0;
57  for (const auto& lcId : lcs) {
58  const auto fraction = 1.f / simTracksters[stId].vertex_multiplicity(lcInSimTst++);
59 
60  const auto lc_find_it = lcToSimTracksterId_Map.find(lcId);
61  if (lc_find_it == lcToSimTracksterId_Map.end()) {
62  lcToSimTracksterId_Map[lcId] = std::vector<hgcal::lcInfoInTrackster>();
63  }
64  lcToSimTracksterId_Map[lcId].emplace_back(stId, fraction);
65 
66  tssInSimTrackster[stId].energy += fraction * layerClusters[lcId].energy();
67  tssInSimTrackster[stId].lcs_and_fractions.emplace_back(lcId, fraction);
68  }
69  } // end of loop over SimTracksters
70 
71 #ifdef EDM_ML_DEBUG
72  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
73  << "tssInSimTrackster INFO (Only SimTrackster filled at the moment)" << std::endl;
74  for (size_t st = 0; st < tssInSimTrackster.size(); ++st) {
75  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "For SimTrackster Idx: " << st << " we have: " << std::endl;
76  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
77  << "\tSimTracksterIdx:\t" << tssInSimTrackster[st].simTracksterId << std::endl;
78  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "\tEnergy:\t" << tssInSimTrackster[st].energy << std::endl;
79  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "\t# of clusters:\t" << layerClusters.size() << std::endl;
80  double tot_energy = 0.;
81  for (auto const& lcaf : tssInSimTrackster[st].lcs_and_fractions) {
82  const auto lcEn = lcaf.second * layerClusters[lcaf.first].energy();
83  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
84  << "\tLC/fraction/energy: " << (uint32_t)lcaf.first << "/" << lcaf.second << "/" << lcEn << std::endl;
85  tot_energy += lcEn;
86  }
87  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "\tTot Sum lcaf: " << tot_energy << std::endl;
88  for (auto const& ts : tssInSimTrackster[st].tracksterIdToEnergyAndScore) {
89  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
90  << "\ttsIdx/energy/score: " << ts.first << "/" << ts.second.first << "/" << ts.second.second << std::endl;
91  }
92  }
93 
94  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "lcToSimTracksterId_Map INFO" << std::endl;
95  for (auto const& lc : lcToSimTracksterId_Map) {
96  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
97  << "For lcId: " << (uint32_t)lc.first
98  << " we have found the following connections with SimTracksters:" << std::endl;
99  for (auto const& st : lc.second) {
100  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
101  << "\tSimTrackster Id: " << st.clusterId << " with fraction: " << st.fraction
102  << " and energy: " << st.fraction * layerClusters[lc.first].energy() << std::endl;
103  }
104  }
105 #endif
106 
107  // Fill stsInTrackster and update tssInSimTrackster
108  // this contains the ids of the simTracksters contributing with at least one
109  // hit to the Trackster. To be returned since this contains the information
110  // to compute the Trackster-To-SimTrackster score.
111  hgcal::tracksterToSimTrackster stsInTrackster; // tsId->(stId,score)
112  stsInTrackster.resize(nTracksters);
113 
114  for (unsigned int tsId = 0; tsId < nTracksters; ++tsId) {
115  for (unsigned int i = 0; i < tracksters[tsId].vertices().size(); ++i) {
116  const auto lcId = tracksters[tsId].vertices(i);
117  const auto lcFractionInTs = 1.f / tracksters[tsId].vertex_multiplicity(i);
118 
119  const auto lc_find_in_ST = lcToSimTracksterId_Map.find(lcId);
120 
121  if (lc_find_in_ST != lcToSimTracksterId_Map.end()) {
122  //Loops through all the simTracksters that contain the layer cluster under study
123  //Here is time to update the tssInSimTrackster and connect the SimTrackster with all
124  //the Tracksters that have the current layer cluster matched.
125  for (const auto& st : lc_find_in_ST->second) {
126  //tssInSimTrackster[simTracksterId][layerclusterId]-> (energy,score)
127  //ST_i - > TS_j, TS_k, ...
128  tssInSimTrackster[st.clusterId].tracksterIdToEnergyAndScore[tsId].first +=
129  lcFractionInTs * st.fraction * layerClusters[lcId].energy();
130  //TS_i -> ST_j, ST_k, ...
131  stsInTrackster[tsId].emplace_back(st.clusterId, std::make_pair(0.f, 0.f));
132  }
133  }
134  } // End loop over LayerClusters in Trackster
135  } // End of loop over Tracksters
136 
137 #ifdef EDM_ML_DEBUG
138  for (unsigned int tsId = 0; tsId < nTracksters; ++tsId) {
139  for (const auto& lcId : tracksters[tsId].vertices()) {
140  const auto& hits_and_fractions = layerClusters[lcId].hitsAndFractions();
141  unsigned int numberOfHitsInLC = hits_and_fractions.size();
142 
143  // This vector will store, for each hit in the Layercluster, the index of
144  // the SimTrackster that contributed the most, in terms of energy, to it.
145  // Special values are:
146  //
147  // -2 --> the reconstruction fraction of the RecHit is 0 (used in the past to monitor Halo Hits)
148  // -3 --> same as before with the added condition that no SimTrackster has been linked to this RecHit
149  // -1 --> the reco fraction is >0, but no SimTrackster has been linked to it
150  // >=0 --> index of the linked SimTrackster
151  std::vector<int> hitsToSimTracksterId(numberOfHitsInLC);
152  // This will store the index of the SimTrackster linked to the LayerCluster that has the largest number of hits in common.
153  int maxSTId_byNumberOfHits = -1;
154  // This will store the maximum number of shared hits between a LayerCluster and a SimTrackster
155  unsigned int maxSTNumberOfHitsInLC = 0;
156  // This will store the index of the SimTrackster linked to the LayerCluster that has the largest energy in common.
157  int maxSTId_byEnergy = -1;
158  // This will store the maximum number of shared energy between a LayerCluster and a SimTrackster
159  float maxEnergySharedLCandST = 0.f;
160  // This will store the fraction of the LayerCluster energy shared with the best(energy) SimTrackster: e_shared/lc_energy
161  float energyFractionOfLCinST = 0.f;
162  // This will store the fraction of the SimTrackster energy shared with the Trackster: e_shared/sc_energy
163  float energyFractionOfSTinLC = 0.f;
164  std::unordered_map<unsigned, unsigned> occurrencesSTinLC;
165  unsigned int numberOfNoiseHitsInLC = 0;
166  std::unordered_map<unsigned, float> STEnergyInLC;
167 
168  const auto lc_find_in_ST = lcToSimTracksterId_Map.find(lcId);
169  for (unsigned int hitId = 0; hitId < numberOfHitsInLC; hitId++) {
170  const auto rhFraction = hits_and_fractions[hitId].second;
171  // if the fraction is zero or the hit does not belong to any SimTrackster,
172  // set the SimTrackster Id for the hit to -1; this will
173  // contribute to the number of noise hits
174 
175  // MR Remove the case in which the fraction is 0, since this could be a
176  // real hit that has been marked as halo.
177  if (rhFraction == 0.) {
178  hitsToSimTracksterId[hitId] = -2;
179  }
180  //Now check if there are SimTracksters linked to this rechit of the layercluster
181  if (lc_find_in_ST == lcToSimTracksterId_Map.end()) {
182  hitsToSimTracksterId[hitId] -= 1;
183  } else {
184  auto maxSTEnergyInLC = 0.f;
185  auto maxSTId = -1;
186  //Loop through all the linked SimTracksters
187  for (const auto& st : lc_find_in_ST->second) {
188  const auto stId = st.clusterId;
189  STEnergyInLC[stId] += st.fraction * layerClusters[lcId].energy();
190  // Keep track of which SimTrackster contributed the most, in terms
191  // of energy, to this specific Layer Cluster.
192  if (STEnergyInLC[stId] > maxSTEnergyInLC) {
193  maxSTEnergyInLC = STEnergyInLC[stId];
194  maxSTId = stId;
195  }
196  }
197  hitsToSimTracksterId[hitId] = maxSTId;
198  }
199  } // End loop over hits on a LayerCluster
200 
201  for (const auto& c : hitsToSimTracksterId) {
202  if (c < 0) {
203  numberOfNoiseHitsInLC++;
204  } else {
205  occurrencesSTinLC[c]++;
206  }
207  }
208 
209  for (const auto& c : occurrencesSTinLC) {
210  if (c.second > maxSTNumberOfHitsInLC) {
211  maxSTId_byNumberOfHits = c.first;
212  maxSTNumberOfHitsInLC = c.second;
213  }
214  }
215 
216  for (const auto& c : STEnergyInLC) {
217  if (c.second > maxEnergySharedLCandST) {
218  maxSTId_byEnergy = c.first;
219  maxEnergySharedLCandST = c.second;
220  }
221  }
222 
223  float totalSTEnergyOnLayer = 0.f;
224  if (maxSTId_byEnergy >= 0) {
225  totalSTEnergyOnLayer = tssInSimTrackster[maxSTId_byEnergy].energy;
226  energyFractionOfSTinLC = maxEnergySharedLCandST / totalSTEnergyOnLayer;
227  if (tracksters[tsId].raw_energy() > 0.f) {
228  energyFractionOfLCinST = maxEnergySharedLCandST / tracksters[tsId].raw_energy();
229  }
230  }
231 
232  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
233  << std::setw(12) << "TracksterID:\t" << std::setw(12) << "layerCluster\t" << std::setw(10) << "lc energy\t"
234  << std::setw(5) << "nhits\t" << std::setw(12) << "noise hits\t" << std::setw(22) << "maxSTId_byNumberOfHits\t"
235  << std::setw(8) << "nhitsST\t" << std::setw(13) << "maxSTId_byEnergy\t" << std::setw(20)
236  << "maxEnergySharedLCandST\t" << std::setw(22) << "totalSTEnergyOnLayer\t" << std::setw(22)
237  << "energyFractionOfLCinST\t" << std::setw(25) << "energyFractionOfSTinLC\t"
238  << "\n";
239  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
240  << std::setw(12) << tsId << "\t" << std::setw(12) << lcId << "\t" << std::setw(10)
241  << tracksters[tsId].raw_energy() << "\t" << std::setw(5) << numberOfHitsInLC << "\t" << std::setw(12)
242  << numberOfNoiseHitsInLC << "\t" << std::setw(22) << maxSTId_byNumberOfHits << "\t" << std::setw(8)
243  << maxSTNumberOfHitsInLC << "\t" << std::setw(13) << maxSTId_byEnergy << "\t" << std::setw(20)
244  << maxEnergySharedLCandST << "\t" << std::setw(22) << totalSTEnergyOnLayer << "\t" << std::setw(22)
245  << energyFractionOfLCinST << "\t" << std::setw(25) << energyFractionOfSTinLC << "\n";
246  } // End of loop over LayerClusters in Trackster
247  } // End of loop over Tracksters
248 
249  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
250  << "Improved tssInSimTrackster INFO (Now containing the linked tracksters id and energy - score still empty)"
251  << std::endl;
252  for (size_t st = 0; st < tssInSimTrackster.size(); ++st) {
253  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "For SimTrackster Idx: " << st << " we have: " << std::endl;
254  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
255  << " SimTracksterIdx: " << tssInSimTrackster[st].simTracksterId << std::endl;
256  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "\tEnergy:\t" << tssInSimTrackster[st].energy << std::endl;
257  double tot_energy = 0.;
258  for (auto const& lcaf : tssInSimTrackster[st].lcs_and_fractions) {
259  const auto lcEn = lcaf.second * layerClusters[lcaf.first].energy();
260  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
261  << "\tLC/fraction/energy: " << (uint32_t)lcaf.first << "/" << lcaf.second << "/" << lcEn << std::endl;
262  tot_energy += lcEn;
263  }
264  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "\tTot Sum lcaf: " << tot_energy << std::endl;
265  for (auto const& ts : tssInSimTrackster[st].tracksterIdToEnergyAndScore) {
266  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
267  << "\ttsIdx/energy/score: " << ts.first << "/" << ts.second.first << "/" << ts.second.second << std::endl;
268  }
269  }
270 
271  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "Improved lcToSimTracksterId_Map INFO" << std::endl;
272  for (auto const& lc : lcToSimTracksterId_Map) {
273  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
274  << "For lcId: " << (uint32_t)lc.first
275  << " we have found the following connections with SimTracksters:" << std::endl;
276  for (auto const& st : lc.second) {
277  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
278  << " SimTrackster Id: " << st.clusterId << " with fraction: " << st.fraction
279  << " and energy: " << st.fraction * layerClusters[lc.first].energy() << std::endl;
280  }
281  }
282 #endif
283 
284  // Update stsInTrackster; compute the score Trackster-to-SimTrackster,
285  // together with the returned AssociationMap
286  for (unsigned int tsId = 0; tsId < nTracksters; ++tsId) {
287  // The SimTracksters contributing to the Trackster's LayerClusters should already be unique.
288  // find the unique SimTracksters id contributing to the Trackster's LayerClusters
289  std::sort(stsInTrackster[tsId].begin(), stsInTrackster[tsId].end());
290  auto last = std::unique(stsInTrackster[tsId].begin(), stsInTrackster[tsId].end());
291  stsInTrackster[tsId].erase(last, stsInTrackster[tsId].end());
292 
293  // If a reconstructed Trackster has energy 0 but is linked to a
294  // SimTrackster, assigned score 1
295  if (tracksters[tsId].raw_energy() == 0. && !stsInTrackster[tsId].empty()) {
296  for (auto& stPair : stsInTrackster[tsId]) {
297  stPair.second.second = 1.;
298  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
299  << "TracksterId:\t " << tsId << "\tST id:\t" << stPair.first << "\tenergy" << stPair.second.first
300  << "\tscore\t " << stPair.second.second << "\n";
301  }
302  continue;
303  }
304 
305  float invTracksterEnergyWeight = 0.f;
306  for (unsigned int i = 0; i < tracksters[tsId].vertices().size(); ++i) {
307  const auto lcId = tracksters[tsId].vertices(i);
308  const auto lcFractionInTs = 1.f / tracksters[tsId].vertex_multiplicity(i);
309 
310  const auto& hits_and_fractions = layerClusters[lcId].hitsAndFractions();
311  // Compute the correct normalization
312  for (auto const& haf : hits_and_fractions) {
313  invTracksterEnergyWeight += std::pow(lcFractionInTs * haf.second * hitMap_->at(haf.first)->energy(), 2);
314  }
315  }
316  invTracksterEnergyWeight = 1.f / invTracksterEnergyWeight;
317 
318  for (unsigned int i = 0; i < tracksters[tsId].vertices().size(); ++i) {
319  const auto lcId = tracksters[tsId].vertices(i);
320  const auto lcFractionInTs = 1.f / tracksters[tsId].vertex_multiplicity(i);
321 
322  const bool lcWithST = (lcToSimTracksterId_Map.find(lcId) != lcToSimTracksterId_Map.end());
323 
324  float lcEnergyWeight = pow(layerClusters[lcId].energy(), 2);
325 
326  for (auto& stPair : stsInTrackster[tsId]) {
327  float stFraction = 0.f;
328  if (lcWithST) {
329  const auto findLCIt = std::find(lcToSimTracksterId_Map[lcId].begin(),
330  lcToSimTracksterId_Map[lcId].end(),
331  hgcal::lcInfoInTrackster{stPair.first, 0.f});
332  if (findLCIt != lcToSimTracksterId_Map[lcId].end()) {
333  stFraction = findLCIt->fraction;
334  }
335  }
336  stPair.second.second +=
337  (lcFractionInTs - stFraction) * (lcFractionInTs - stFraction) * lcEnergyWeight * invTracksterEnergyWeight;
338 #ifdef EDM_ML_DEBUG
339  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
340  << "lcId:\t" << (uint32_t)lcId << "\ttracksterId:\t" << tsId << "\ttsFraction,stFraction:\t"
341  << lcFractionInTs << ", " << stFraction << "\tlcEnergyWeight:\t" << lcEnergyWeight << "\tcurrent score:\t"
342  << stPair.second.second << "\tinvTracksterEnergyWeight:\t" << invTracksterEnergyWeight << "\n";
343 #endif
344  }
345  } // End of loop over LayerClusters in Trackster
346 
347 #ifdef EDM_ML_DEBUG
348  if (stsInTrackster[tsId].empty())
349  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "trackster Id:\t" << tsId << "\tST id:\t-1"
350  << "\tscore\t-1\n";
351 #endif
352  } // End of loop over Tracksters
353 
354  // Compute the SimTrackster-To-Trackster score
355  for (const auto& stId : sTIndices) {
356  float invSTEnergyWeight = 0.f;
357 
358  const unsigned int STNumberOfLCs = tssInSimTrackster[stId].lcs_and_fractions.size();
359  if (STNumberOfLCs == 0)
360  continue;
361 
362 #ifdef EDM_ML_DEBUG
363  int tsWithMaxEnergyInST = -1;
364  //energy of the most energetic TS from all that were linked to ST
365  float maxEnergyTSinST = 0.f;
366  float STenergy = tssInSimTrackster[stId].energy;
367  //most energetic TS from all TSs linked to ST over ST energy.
368  float STEnergyFractionInTS = 0.f;
369  for (const auto& ts : tssInSimTrackster[stId].tracksterIdToEnergyAndScore) {
370  if (ts.second.first > maxEnergyTSinST) {
371  maxEnergyTSinST = ts.second.first;
372  tsWithMaxEnergyInST = ts.first;
373  }
374  }
375  if (STenergy > 0.f)
376  STEnergyFractionInTS = maxEnergyTSinST / STenergy;
377 
378  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
379  << std::setw(12) << "simTrackster\t" << std::setw(15) << "st total energy\t" << std::setw(15)
380  << "stEnergyOnLayer\t" << std::setw(14) << "STNhitsOnLayer\t" << std::setw(18) << "tsWithMaxEnergyInST\t"
381  << std::setw(15) << "maxEnergyTSinST\t" << std::setw(20) << "STEnergyFractionInTS"
382  << "\n";
383  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
384  << std::setw(12) << stId << "\t" << std::setw(15) << simTracksters[stId].raw_energy() << "\t" << std::setw(15)
385  << STenergy << "\t" << std::setw(14) << STNumberOfLCs << "\t" << std::setw(18) << tsWithMaxEnergyInST << "\t"
386  << std::setw(15) << maxEnergyTSinST << "\t" << std::setw(20) << STEnergyFractionInTS << "\n";
387 #endif
388 
389  // Compute the correct normalization
390  for (auto const& lcaf : tssInSimTrackster[stId].lcs_and_fractions) {
391  invSTEnergyWeight += std::pow(lcaf.second * layerClusters[lcaf.first].energy(), 2);
392  }
393  invSTEnergyWeight = 1.f / invSTEnergyWeight;
394 
395  for (unsigned int i = 0; i < STNumberOfLCs; ++i) {
396  auto& st_lcId = tssInSimTrackster[stId].lcs_and_fractions[i].first;
397  auto& st_lcFraction = tssInSimTrackster[stId].lcs_and_fractions[i].second;
398 
399  if (st_lcFraction == 0.f)
400  continue; // hopefully this should never happen
401  float lcEnergyWeight = pow(layerClusters[st_lcId].energy(), 2);
402  for (auto& tsPair : tssInSimTrackster[stId].tracksterIdToEnergyAndScore) {
403  unsigned int tsId = tsPair.first;
404 
405  for (unsigned int i = 0; i < tracksters[tsId].vertices().size(); ++i) {
406  const auto tsFraction = 1.f / tracksters[tsId].vertex_multiplicity(i);
407 
408  tsPair.second.second +=
409  (tsFraction - st_lcFraction) * (tsFraction - st_lcFraction) * lcEnergyWeight * invSTEnergyWeight;
410 #ifdef EDM_ML_DEBUG
411  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
412  << "STLCId:\t" << (uint32_t)st_lcId << "\tTracksterId:\t" << tsId << "\t"
413  << "tsFraction, stFraction:\t" << tsFraction << ", " << st_lcFraction << "\t"
414  << "lcEnergyWeight:\t" << lcEnergyWeight << "\t"
415  << "current score:\t" << tsPair.second.second << "\t"
416  << "invSTEnergyWeight:\t" << invSTEnergyWeight << "\n";
417 #endif
418  } // End of loop over Trackster's LayerClusters
419  } // End of loop over Tracksters linked to hits of this SimTrackster
420  } // End of loop over hits of SimTrackster on a Layer
421 #ifdef EDM_ML_DEBUG
422  if (tssInSimTrackster[stId].tracksterIdToEnergyAndScore.empty())
423  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "ST Id:\t" << stId << "\tTS id:\t-1 "
424  << "\tscore\t-1\n";
425 
426  for (const auto& tsPair : tssInSimTrackster[stId].tracksterIdToEnergyAndScore) {
427  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl")
428  << "ST Id: \t" << stId << "\t TS id: \t" << tsPair.first << "\t score \t" << tsPair.second.second
429  << "\t shared energy:\t" << tsPair.second.first << "\t shared energy fraction:\t"
430  << (tsPair.second.first / STenergy) << "\n";
431  }
432 #endif
433  } // End loop over SimTrackster indices
434  return {stsInTrackster, tssInSimTrackster};
435 }
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
def unique(seq, keepstr=True)
Definition: tier0.py:24
double f[11][100]
std::vector< hgcal::simTracksterOnLayer > simTracksterToTrackster
std::vector< std::vector< std::pair< unsigned int, std::pair< float, float > > > > tracksterToSimTrackster
const std::unordered_map< DetId, const HGCRecHit * > * hitMap_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
#define LogDebug(id)

Member Data Documentation

◆ hardScatterOnly_

const bool TSToSimTSAssociatorByEnergyScoreImpl::hardScatterOnly_
private

Definition at line 84 of file TSToSimTSAssociatorByEnergyScoreImpl.h.

◆ hitMap_

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

Definition at line 86 of file TSToSimTSAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ layers_

unsigned TSToSimTSAssociatorByEnergyScoreImpl::layers_
private

◆ productGetter_

edm::EDProductGetter const* TSToSimTSAssociatorByEnergyScoreImpl::productGetter_
private

◆ recHitTools_

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