CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
LCToSCAssociatorByEnergyScoreImpl< HIT > Class Template Reference

#include <LCToSCAssociatorByEnergyScoreImpl.h>

Inheritance diagram for LCToSCAssociatorByEnergyScoreImpl< HIT >:
ticl::LayerClusterToSimClusterAssociatorBaseImpl

Public Member Functions

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

Private Member Functions

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

Private Attributes

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

Detailed Description

template<typename HIT>
class LCToSCAssociatorByEnergyScoreImpl< HIT >

Definition at line 65 of file LCToSCAssociatorByEnergyScoreImpl.h.

Constructor & Destructor Documentation

◆ LCToSCAssociatorByEnergyScoreImpl()

template<typename HIT >
LCToSCAssociatorByEnergyScoreImpl< HIT >::LCToSCAssociatorByEnergyScoreImpl ( edm::EDProductGetter const &  productGetter,
bool  hardScatterOnly,
std::shared_ptr< hgcal::RecHitTools recHitTools,
const std::unordered_map< DetId, const unsigned int > *  hitMap,
std::vector< const HIT *> &  hits 
)
explicit

Definition at line 7 of file LCToSCAssociatorByEnergyScoreImpl.cc.

References ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), LCToSCAssociatorByEnergyScoreImpl< HIT >::layers_, and LCToSCAssociatorByEnergyScoreImpl< HIT >::recHitTools_.

13  : hardScatterOnly_(hardScatterOnly),
14  recHitTools_(recHitTools),
15  hitMap_(hitMap),
17  hits_(hits) {
18  if constexpr (std::is_same_v<HIT, HGCRecHit>)
19  layers_ = recHitTools_->lastLayerBH();
20  else
21  layers_ = 6; //EB + 4 HB + HO
22 }
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()

template<typename HIT >
ticl::RecoToSimCollectionWithSimClusters LCToSCAssociatorByEnergyScoreImpl< HIT >::associateRecoToSim ( const edm::Handle< reco::CaloClusterCollection > &  cCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
overridevirtual

Associate a LayerCluster to SimClusters.

Reimplemented from ticl::LayerClusterToSimClusterAssociatorBaseImpl.

Definition at line 516 of file LCToSCAssociatorByEnergyScoreImpl.cc.

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

517  {
519  const auto& links = makeConnections(cCCH, sCCH);
520 
521  const auto& scsInLayerCluster = std::get<0>(links);
522  for (size_t lcId = 0; lcId < scsInLayerCluster.size(); ++lcId) {
523  for (auto& scPair : scsInLayerCluster[lcId]) {
524  LogDebug("LCToSCAssociatorByEnergyScoreImpl")
525  << "layerCluster Id: \t" << lcId << "\t SC id: \t" << scPair.first << "\t score \t" << scPair.second << "\n";
526  // Fill AssociationMap
527  returnValue.insert(edm::Ref<reco::CaloClusterCollection>(cCCH, lcId), // Ref to LC
528  std::make_pair(edm::Ref<SimClusterCollection>(sCCH, scPair.first),
529  scPair.second) // Pair <Ref to SC, score>
530  );
531  }
532  }
533  return returnValue;
534 }
#define LogDebug(id)
ticl::association makeConnections(const edm::Handle< reco::CaloClusterCollection > &cCH, const edm::Handle< SimClusterCollection > &sCCH) const

◆ associateSimToReco()

template<typename HIT >
ticl::SimToRecoCollectionWithSimClusters LCToSCAssociatorByEnergyScoreImpl< HIT >::associateSimToReco ( const edm::Handle< reco::CaloClusterCollection > &  cCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
overridevirtual

Associate a SimCluster to LayerClusters.

Reimplemented from ticl::LayerClusterToSimClusterAssociatorBaseImpl.

Definition at line 537 of file LCToSCAssociatorByEnergyScoreImpl.cc.

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

538  {
540  const auto& links = makeConnections(cCCH, sCCH);
541  const auto& lcsInSimCluster = std::get<1>(links);
542  for (size_t scId = 0; scId < lcsInSimCluster.size(); ++scId) {
543  for (size_t layerId = 0; layerId < lcsInSimCluster[scId].size(); ++layerId) {
544  for (auto& lcPair : lcsInSimCluster[scId][layerId].layerClusterIdToEnergyAndScore) {
545  returnValue.insert(
546  edm::Ref<SimClusterCollection>(sCCH, scId), // Ref to SC
547  std::make_pair(edm::Ref<reco::CaloClusterCollection>(cCCH, lcPair.first), // Pair <Ref to LC,
548  std::make_pair(lcPair.second.first, lcPair.second.second)) // pair <energy, score> >
549  );
550  }
551  }
552  }
553  return returnValue;
554 }
ticl::association makeConnections(const edm::Handle< reco::CaloClusterCollection > &cCH, const edm::Handle< SimClusterCollection > &sCCH) const

◆ makeConnections()

template<typename HIT >
ticl::association LCToSCAssociatorByEnergyScoreImpl< HIT >::makeConnections ( const edm::Handle< reco::CaloClusterCollection > &  cCH,
const edm::Handle< SimClusterCollection > &  sCCH 
) const
private

Definition at line 25 of file LCToSCAssociatorByEnergyScoreImpl.cc.

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

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

template<typename HIT >
const bool LCToSCAssociatorByEnergyScoreImpl< HIT >::hardScatterOnly_
private

Definition at line 82 of file LCToSCAssociatorByEnergyScoreImpl.h.

◆ hitMap_

template<typename HIT >
const std::unordered_map<DetId, const unsigned int>* LCToSCAssociatorByEnergyScoreImpl< HIT >::hitMap_
private

Definition at line 84 of file LCToSCAssociatorByEnergyScoreImpl.h.

◆ hits_

template<typename HIT >
std::vector<const HIT *> LCToSCAssociatorByEnergyScoreImpl< HIT >::hits_
private

Definition at line 89 of file LCToSCAssociatorByEnergyScoreImpl.h.

◆ layers_

template<typename HIT >
unsigned LCToSCAssociatorByEnergyScoreImpl< HIT >::layers_
private

◆ productGetter_

template<typename HIT >
edm::EDProductGetter const* LCToSCAssociatorByEnergyScoreImpl< HIT >::productGetter_
private

Definition at line 86 of file LCToSCAssociatorByEnergyScoreImpl.h.

◆ recHitTools_

template<typename HIT >
std::shared_ptr<hgcal::RecHitTools> LCToSCAssociatorByEnergyScoreImpl< HIT >::recHitTools_
private