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

Public Member Functions

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

Private Member Functions

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

Definition at line 5 of file TSToSimTSAssociatorByEnergyScoreImpl.cc.

References layers_, and recHitTools_.

11  : hardScatterOnly_(hardScatterOnly),
12  recHitTools_(recHitTools),
13  hitMap_(hitMap),
14  hits_(hits),
16  layers_ = recHitTools_->lastLayerBH();
17 }
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::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 ticl::TracksterToSimTracksterAssociatorBaseImpl.

Definition at line 443 of file TSToSimTSAssociatorByEnergyScoreImpl.cc.

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

446  {
448  const auto& links = makeConnections(tCH, lCCH, sTCH);
449 
450  const auto& stsInTrackster = std::get<0>(links);
451  for (size_t tsId = 0; tsId < stsInTrackster.size(); ++tsId) {
452  for (auto& stPair : stsInTrackster[tsId]) {
453  LogDebug("TSToSimTSAssociatorByEnergyScoreImpl") << "Trackster Id:\t" << tsId << "\tSimTrackster id:\t"
454  << stPair.first << "\tscore:\t" << stPair.second.second << "\n";
455  // Fill AssociationMap
456  returnValue.insert(
457  edm::Ref<ticl::TracksterCollection>(tCH, tsId), // Ref to TS
458  std::make_pair(edm::Ref<ticl::TracksterCollection>(sTCH, stPair.first), //Pair <Refo to TS>
459  std::make_pair(stPair.second.first, stPair.second.second)) // Pair <energy, score>
460  );
461  }
462  }
463  return returnValue;
464 }
ticl::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()

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

Definition at line 466 of file TSToSimTSAssociatorByEnergyScoreImpl.cc.

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

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

◆ makeConnections()

ticl::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 19 of file TSToSimTSAssociatorByEnergyScoreImpl.cc.

References DummyCfis::c, relativeConstraints::empty, mps_fire::end, hcalRecHitTable_cff::energy, f, spr::find(), cms::cuda::for(), HLT_2024v14_cff::fraction, hitMap_, hits_, 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().

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

◆ hitMap_

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

Definition at line 87 of file TSToSimTSAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ hits_

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

Definition at line 88 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