CMS 3D CMS Logo

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

#include <MultiClusterAssociatorByEnergyScoreImpl.h>

Inheritance diagram for MultiClusterAssociatorByEnergyScoreImpl:
hgcal::MultiClusterToCaloParticleAssociatorBaseImpl

Public Member Functions

hgcal::RecoToSimCollectionWithMultiClusters associateRecoToSim (const edm::Handle< reco::HGCalMultiClusterCollection > &mCCH, const edm::Handle< CaloParticleCollection > &cPCH) const override
 Associate a MultiCluster to CaloParticles. More...
 
hgcal::SimToRecoCollectionWithMultiClusters associateSimToReco (const edm::Handle< reco::HGCalMultiClusterCollection > &mCCH, const edm::Handle< CaloParticleCollection > &cPCH) const override
 Associate a CaloParticle to MultiClusters. More...
 
 MultiClusterAssociatorByEnergyScoreImpl (edm::EDProductGetter const &, bool, std::shared_ptr< hgcal::RecHitTools >, const std::unordered_map< DetId, const HGCRecHit *> *&)
 
- Public Member Functions inherited from hgcal::MultiClusterToCaloParticleAssociatorBaseImpl
 MultiClusterToCaloParticleAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~MultiClusterToCaloParticleAssociatorBaseImpl ()
 Destructor. More...
 

Private Member Functions

hgcal::association makeConnections (const edm::Handle< reco::HGCalMultiClusterCollection > &mCCH, const edm::Handle< CaloParticleCollection > &cPCH) 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 47 of file MultiClusterAssociatorByEnergyScoreImpl.h.

Constructor & Destructor Documentation

◆ MultiClusterAssociatorByEnergyScoreImpl()

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

Definition at line 13 of file MultiClusterAssociatorByEnergyScoreImpl.cc.

References layers_, and recHitTools_.

18  : hardScatterOnly_(hardScatterOnly), recHitTools_(recHitTools), hitMap_(hitMap), productGetter_(&productGetter) {
19  layers_ = recHitTools_->lastLayerBH();
20 }
const std::unordered_map< DetId, const HGCRecHit * > * hitMap_
EDProductGetter const * productGetter(std::atomic< void const *> const &iCache)

Member Function Documentation

◆ associateRecoToSim()

hgcal::RecoToSimCollectionWithMultiClusters MultiClusterAssociatorByEnergyScoreImpl::associateRecoToSim ( const edm::Handle< reco::HGCalMultiClusterCollection > &  cCH,
const edm::Handle< CaloParticleCollection > &  cPCH 
) const
overridevirtual

Associate a MultiCluster to CaloParticles.

Reimplemented from hgcal::MultiClusterToCaloParticleAssociatorBaseImpl.

Definition at line 494 of file MultiClusterAssociatorByEnergyScoreImpl.cc.

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

495  {
497  const auto& links = makeConnections(mCCH, cPCH);
498 
499  const auto& cpsInMultiCluster = std::get<0>(links);
500  for (size_t mcId = 0; mcId < cpsInMultiCluster.size(); ++mcId) {
501  for (auto& cpPair : cpsInMultiCluster[mcId]) {
502  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
503  << "multiCluster Id: \t" << mcId << "\t CP id: \t" << cpPair.first << "\t score \t" << cpPair.second << "\n";
504  // Fill AssociationMap
505  returnValue.insert(edm::Ref<reco::HGCalMultiClusterCollection>(mCCH, mcId), // Ref to MC
506  std::make_pair(edm::Ref<CaloParticleCollection>(cPCH, cpPair.first),
507  cpPair.second) // Pair <Ref to CP, score>
508  );
509  }
510  }
511  return returnValue;
512 }
hgcal::association makeConnections(const edm::Handle< reco::HGCalMultiClusterCollection > &mCCH, const edm::Handle< CaloParticleCollection > &cPCH) const
#define LogDebug(id)

◆ associateSimToReco()

hgcal::SimToRecoCollectionWithMultiClusters MultiClusterAssociatorByEnergyScoreImpl::associateSimToReco ( const edm::Handle< reco::HGCalMultiClusterCollection > &  cCH,
const edm::Handle< CaloParticleCollection > &  cPCH 
) const
overridevirtual

Associate a CaloParticle to MultiClusters.

Reimplemented from hgcal::MultiClusterToCaloParticleAssociatorBaseImpl.

Definition at line 514 of file MultiClusterAssociatorByEnergyScoreImpl.cc.

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

515  {
517  const auto& links = makeConnections(mCCH, cPCH);
518  const auto& cPOnLayer = std::get<1>(links);
519  for (size_t cpId = 0; cpId < cPOnLayer.size(); ++cpId) {
520  for (size_t layerId = 0; layerId < cPOnLayer[cpId].size(); ++layerId) {
521  for (auto& mcPair : cPOnLayer[cpId][layerId].multiClusterIdToEnergyAndScore) {
522  returnValue.insert(
523  edm::Ref<CaloParticleCollection>(cPCH, cpId), // Ref to CP
524  std::make_pair(edm::Ref<reco::HGCalMultiClusterCollection>(mCCH, mcPair.first), // Pair <Ref to MC,
525  std::make_pair(mcPair.second.first, mcPair.second.second)) // pair <energy, score> >
526  );
527  }
528  }
529  }
530  return returnValue;
531 }
hgcal::association makeConnections(const edm::Handle< reco::HGCalMultiClusterCollection > &mCCH, const edm::Handle< CaloParticleCollection > &cPCH) const

◆ makeConnections()

hgcal::association MultiClusterAssociatorByEnergyScoreImpl::makeConnections ( const edm::Handle< reco::HGCalMultiClusterCollection > &  mCCH,
const edm::Handle< CaloParticleCollection > &  cPCH 
) const
private

Definition at line 22 of file MultiClusterAssociatorByEnergyScoreImpl.cc.

References HltBtagPostValidation_cff::c, caloTruthCellsProducer_cfi::caloParticles, CommonMethods::cp(), relativeConstraints::empty, mps_fire::end, HCALHighEnergyHPDFilter_cfi::energy, f, spr::find(), cms::cuda::for(), newFWLiteAna::found, h, hitMap_, SimCluster::hits_and_fractions(), mps_fire::i, dqmiolumiharvest::j, dqmdumpme::last, layers_, LogDebug, CaloTowersParam_cfi::mc, SiStripPI::min, edm::Handle< T >::product(), recHitTools_, removeCPFromPU(), jetUpdater_cfi::sort, tier0::unique(), and findQualityFiles::v.

Referenced by associateRecoToSim(), and associateSimToReco().

23  {
24  // 1. Extract collections and filter CaloParticles, if required
25  const auto& mClusters = *mCCH.product();
26  const auto& caloParticles = *cPCH.product();
27  auto nMultiClusters = mClusters.size();
28  //Consider CaloParticles coming from the hard scatterer, excluding the PU contribution.
29  std::vector<size_t> cPIndices;
30  //Consider CaloParticles coming from the hard scatterer
31  //excluding the PU contribution and save the indices.
32  removeCPFromPU(caloParticles, cPIndices, false);
33  auto nCaloParticles = cPIndices.size();
34 
35  std::vector<size_t> cPSelectedIndices;
36  removeCPFromPU(caloParticles, cPSelectedIndices, true);
37 
38  //cPOnLayer[caloparticle][layer]
39  //This defines a "caloParticle on layer" concept. It is only filled in case
40  //that caloParticle has a reconstructed hit related via detid. So, a cPOnLayer[i][j] connects a
41  //specific caloParticle i in layer j with:
42  //1. the sum of all recHits energy times fraction of the relevant simHit in layer j related to that caloParticle i.
43  //2. the hits and fractions of that caloParticle i in layer j.
44  //3. the layer clusters with matched recHit id.
46  cPOnLayer.resize(nCaloParticles);
47  for (unsigned int i = 0; i < nCaloParticles; ++i) {
48  auto cpIndex = cPIndices[i];
49  cPOnLayer[cpIndex].resize(layers_ * 2);
50  for (unsigned int j = 0; j < layers_ * 2; ++j) {
51  cPOnLayer[cpIndex][j].caloParticleId = cpIndex;
52  cPOnLayer[cpIndex][j].energy = 0.f;
53  cPOnLayer[cpIndex][j].hits_and_fractions.clear();
54  }
55  }
56 
57  std::unordered_map<DetId, std::vector<hgcal::detIdInfoInCluster>> detIdToCaloParticleId_Map;
58  // Fill detIdToCaloParticleId_Map and update cPOnLayer
59  for (const auto& cpId : cPIndices) {
60  //take sim clusters
61  const SimClusterRefVector& simClusterRefVector = caloParticles[cpId].simClusters();
62  //loop through sim clusters
63  for (const auto& it_sc : simClusterRefVector) {
64  const SimCluster& simCluster = (*(it_sc));
65  const auto& hits_and_fractions = simCluster.hits_and_fractions();
66  for (const auto& it_haf : hits_and_fractions) {
67  const auto hitid = (it_haf.first);
68  const auto cpLayerId =
69  recHitTools_->getLayerWithOffset(hitid) + layers_ * ((recHitTools_->zside(hitid) + 1) >> 1) - 1;
70  const auto itcheck = hitMap_->find(hitid);
71  if (itcheck != hitMap_->end()) {
72  //Since the current hit from sim cluster has a reconstructed hit with the same detid,
73  //make a map that will connect a detid with:
74  //1. the caloParticles that have a simcluster with sim hits in that cell via caloParticle id.
75  //2. the sum of all simHits fractions that contributes to that detid.
76  //So, keep in mind that in case of multiple caloParticles contributing in the same cell
77  //the fraction is the sum over all caloParticles. So, something like:
78  //detid: (caloParticle 1, sum of hits fractions in that detid over all cp) , (caloParticle 2, sum of hits fractions in that detid over all cp), (caloParticle 3, sum of hits fractions in that detid over all cp) ...
79  auto hit_find_it = detIdToCaloParticleId_Map.find(hitid);
80  if (hit_find_it == detIdToCaloParticleId_Map.end()) {
81  detIdToCaloParticleId_Map[hitid] = std::vector<hgcal::detIdInfoInCluster>();
82  detIdToCaloParticleId_Map[hitid].emplace_back(cpId, it_haf.second);
83  } else {
84  auto findHitIt = std::find(detIdToCaloParticleId_Map[hitid].begin(),
85  detIdToCaloParticleId_Map[hitid].end(),
86  hgcal::detIdInfoInCluster{cpId, it_haf.second});
87  if (findHitIt != detIdToCaloParticleId_Map[hitid].end()) {
88  findHitIt->fraction += it_haf.second;
89  } else {
90  detIdToCaloParticleId_Map[hitid].emplace_back(cpId, it_haf.second);
91  }
92  }
93  const HGCRecHit* hit = itcheck->second;
94  //Since the current hit from sim cluster has a reconstructed hit with the same detid,
95  //fill the cPOnLayer[caloparticle][layer] object with energy (sum of all recHits energy times fraction
96  //of the relevant simHit) and keep the hit (detid and fraction) that contributed.
97  cPOnLayer[cpId][cpLayerId].energy += it_haf.second * hit->energy();
98  // We need to compress the hits and fractions in order to have a
99  // reasonable score between CP and LC. Imagine, for example, that a
100  // CP has detID X used by 2 SimClusters with different fractions. If
101  // a single LC uses X with fraction 1 and is compared to the 2
102  // contributions separately, it will be assigned a score != 0, which
103  // is wrong.
104  auto& haf = cPOnLayer[cpId][cpLayerId].hits_and_fractions;
105  auto found = std::find_if(
106  std::begin(haf), std::end(haf), [&hitid](const std::pair<DetId, float>& v) { return v.first == hitid; });
107  if (found != haf.end()) {
108  found->second += it_haf.second;
109  } else {
110  cPOnLayer[cpId][cpLayerId].hits_and_fractions.emplace_back(hitid, it_haf.second);
111  }
112  }
113  } // end of loop through simHits
114  } // end of loop through simclusters
115  } // end of loop through caloParticles
116 
117 #ifdef EDM_ML_DEBUG
118  LogDebug("MultiClusterAssociatorByEnergyScoreImpl") << "cPOnLayer INFO" << std::endl;
119  for (size_t cp = 0; cp < cPOnLayer.size(); ++cp) {
120  LogDebug("MultiClusterAssociatorByEnergyScoreImpl") << "For CaloParticle Idx: " << cp << " we have: " << std::endl;
121  for (size_t cpp = 0; cpp < cPOnLayer[cp].size(); ++cpp) {
122  LogDebug("MultiClusterAssociatorByEnergyScoreImpl") << " On Layer: " << cpp << " we have:" << std::endl;
123  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
124  << " CaloParticleIdx: " << cPOnLayer[cp][cpp].caloParticleId << std::endl;
125  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
126  << " Energy: " << cPOnLayer[cp][cpp].energy << std::endl;
127  double tot_energy = 0.;
128  for (auto const& haf : cPOnLayer[cp][cpp].hits_and_fractions) {
129  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
130  << " Hits/fraction/energy: " << (uint32_t)haf.first << "/" << haf.second << "/"
131  << haf.second * hitMap_->at(haf.first)->energy() << std::endl;
132  tot_energy += haf.second * hitMap_->at(haf.first)->energy();
133  }
134  LogDebug("MultiClusterAssociatorByEnergyScoreImpl") << " Tot Sum haf: " << tot_energy << std::endl;
135  for (auto const& mc : cPOnLayer[cp][cpp].multiClusterIdToEnergyAndScore) {
136  LogDebug("MultiClusterAssociatorByEnergyScoreImpl") << " mcIdx/energy/score: " << mc.first << "/"
137  << mc.second.first << "/" << mc.second.second << std::endl;
138  }
139  }
140  }
141 
142  LogDebug("MultiClusterAssociatorByEnergyScoreImpl") << "detIdToCaloParticleId_Map INFO" << std::endl;
143  for (auto const& cp : detIdToCaloParticleId_Map) {
144  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
145  << "For detId: " << (uint32_t)cp.first
146  << " we have found the following connections with CaloParticles:" << std::endl;
147  for (auto const& cpp : cp.second) {
148  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
149  << " CaloParticle Id: " << cpp.clusterId << " with fraction: " << cpp.fraction
150  << " and energy: " << cpp.fraction * hitMap_->at(cp.first)->energy() << std::endl;
151  }
152  }
153 #endif
154 
155  // Fill detIdToMultiClusterId_Map and cpsInMultiCluster; update cPOnLayer
156  std::unordered_map<DetId, std::vector<hgcal::detIdInfoInMultiCluster>> detIdToMultiClusterId_Map;
157 
158  // this contains the ids of the caloParticles contributing with at least one hit to the multiCluster and the reconstruction error
159  //cpsInMultiCluster[multicluster][CPids]
160  //Connects a multiCluster with all related caloParticles.
161  hgcal::multiClusterToCaloParticle cpsInMultiCluster;
162  cpsInMultiCluster.resize(nMultiClusters);
163 
164  //Loop through multiClusters
165  for (unsigned int mcId = 0; mcId < nMultiClusters; ++mcId) {
166  const auto& hits_and_fractions = mClusters[mcId].hitsAndFractions();
167  if (!hits_and_fractions.empty()) {
168  std::unordered_map<unsigned, float> CPEnergyInMCL;
169  int maxCPId_byNumberOfHits = -1;
170  unsigned int maxCPNumberOfHitsInMCL = 0;
171  int maxCPId_byEnergy = -1;
172  float maxEnergySharedMCLandCP = 0.f;
173  float energyFractionOfMCLinCP = 0.f;
174  float energyFractionOfCPinMCL = 0.f;
175 
176  //In case of matched rechit-simhit, so matched
177  //caloparticle-layercluster-multicluster, we count and save the number of
178  //recHits related to the maximum energy CaloParticle out of all
179  //CaloParticles related to that layer cluster and multiCluster.
180 
181  std::unordered_map<unsigned, unsigned> occurrencesCPinMCL;
182  unsigned int numberOfNoiseHitsInMCL = 0;
183  unsigned int numberOfHaloHitsInMCL = 0;
184  unsigned int numberOfHitsInMCL = 0;
185 
186  //number of hits related to that cluster
187  unsigned int numberOfHitsInLC = hits_and_fractions.size();
188  numberOfHitsInMCL += numberOfHitsInLC;
189  std::unordered_map<unsigned, float> CPEnergyInLC;
190 
191  //hitsToCaloParticleId is a vector of ints, one for each recHit of the
192  //layer cluster under study. If negative, there is no simHit from any CaloParticle related.
193  //If positive, at least one CaloParticle has been found with matched simHit.
194  //In more detail:
195  // 1. hitsToCaloParticleId[hitId] = -3
196  // TN: These represent Halo Cells(N) that have not been
197  // assigned to any CaloParticle (hence the T).
198  // 2. hitsToCaloParticleId[hitId] = -2
199  // FN: There represent Halo Cells(N) that have been assigned
200  // to a CaloParticle (hence the F, since those should have not been marked as halo)
201  // 3. hitsToCaloParticleId[hitId] = -1
202  // FP: These represent Real Cells(P) that have not been
203  // assigned to any CaloParticle (hence the F, since these are fakes)
204  // 4. hitsToCaloParticleId[hitId] >= 0
205  // TP There represent Real Cells(P) that have been assigned
206  // to a CaloParticle (hence the T)
207 
208  std::vector<int> hitsToCaloParticleId(numberOfHitsInLC);
209  //det id of the first hit just to make the lcLayerId variable
210  //which maps the layers in -z: 0->51 and in +z: 52->103
211  const auto firstHitDetId = hits_and_fractions[0].first;
212  int lcLayerId = recHitTools_->getLayerWithOffset(firstHitDetId) +
213  layers_ * ((recHitTools_->zside(firstHitDetId) + 1) >> 1) - 1;
214 
215  //Loop though the hits of the layer cluster under study
216  for (unsigned int hitId = 0; hitId < numberOfHitsInLC; hitId++) {
217  const auto rh_detid = hits_and_fractions[hitId].first;
218  const auto rhFraction = hits_and_fractions[hitId].second;
219 
220  //Since the hit is belonging to the layer cluster, it must also be in the recHits map.
221  const auto itcheck = hitMap_->find(rh_detid);
222  const auto hit = itcheck->second;
223 
224  //Make a map that will connect a detid (that belongs to a recHit of the layer cluster under study,
225  //no need to save others) with:
226  //1. the layer clusters that have recHits in that detid
227  //2. the fraction of the recHit of each layer cluster that contributes to that detid.
228  //So, something like:
229  //detid: (layer cluster 1, hit fraction) , (layer cluster 2, hit fraction), (layer cluster 3, hit fraction) ...
230  //here comparing with the caloParticle map above
231  auto hit_find_in_LC = detIdToMultiClusterId_Map.find(rh_detid);
232  if (hit_find_in_LC == detIdToMultiClusterId_Map.end()) {
233  detIdToMultiClusterId_Map[rh_detid] = std::vector<hgcal::detIdInfoInMultiCluster>();
234  }
235  detIdToMultiClusterId_Map[rh_detid].emplace_back(hgcal::detIdInfoInMultiCluster{mcId, mcId, rhFraction});
236 
237  // Check whether the recHit of the layer cluster under study has a sim hit in the same cell
238  auto hit_find_in_CP = detIdToCaloParticleId_Map.find(rh_detid);
239 
240  // If the fraction is zero or the hit does not belong to any calo
241  // particle, set the caloParticleId for the hit to -1 and this will
242  // contribute to the number of noise hits
243  if (rhFraction == 0.) { // this could be a real hit that has been marked as halo
244  hitsToCaloParticleId[hitId] = -2;
245  numberOfHaloHitsInMCL++;
246  }
247  if (hit_find_in_CP == detIdToCaloParticleId_Map.end()) {
248  hitsToCaloParticleId[hitId] -= 1;
249  } else {
250  auto maxCPEnergyInLC = 0.f;
251  auto maxCPId = -1;
252  for (auto& h : hit_find_in_CP->second) {
253  auto shared_fraction = std::min(rhFraction, h.fraction);
254  //We are in the case where there are caloParticles with simHits connected via detid with the recHit under study
255  //So, from all layers clusters, find the recHits that are connected with a caloParticle and save/calculate the
256  //energy of that caloParticle as the sum over all recHits of the recHits energy weighted
257  //by the caloParticle's fraction related to that recHit.
258  CPEnergyInMCL[h.clusterId] += shared_fraction * hit->energy();
259  //Same but for layer clusters for the cell association per layer
260  CPEnergyInLC[h.clusterId] += shared_fraction * hit->energy();
261  //Here cPOnLayer[caloparticle][layer] described above is set
262  //Here for multiClusters with matched recHit, the CP fraction times hit energy is added and saved
263  cPOnLayer[h.clusterId][lcLayerId].multiClusterIdToEnergyAndScore[mcId].first +=
264  shared_fraction * hit->energy();
265  cPOnLayer[h.clusterId][lcLayerId].multiClusterIdToEnergyAndScore[mcId].second = FLT_MAX;
266  //cpsInMultiCluster[multicluster][CPids]
267  //Connects a multiCluster with all related caloParticles
268  cpsInMultiCluster[mcId].emplace_back(h.clusterId, FLT_MAX);
269  //From all CaloParticles related to a layer cluster, we save id and energy of the caloParticle
270  //that after simhit-rechit matching in layer has the maximum energy.
271  if (shared_fraction > maxCPEnergyInLC) {
272  //energy is used only here. cpid is saved for multiClusters
273  maxCPEnergyInLC = CPEnergyInLC[h.clusterId];
274  maxCPId = h.clusterId;
275  }
276  }
277  //Keep in mind here maxCPId could be zero. So, below ask for negative not including zero to count noise.
278  hitsToCaloParticleId[hitId] = maxCPId;
279  }
280 
281  } //end of loop through recHits of the layer cluster.
282 
283  //Loop through all recHits to count how many of them are noise and how many are matched.
284  //In case of matched rechit-simhit, we count and save the number of recHits related to the maximum energy CaloParticle.
285  for (auto c : hitsToCaloParticleId) {
286  if (c < 0) {
287  numberOfNoiseHitsInMCL++;
288  } else {
289  occurrencesCPinMCL[c]++;
290  }
291  }
292 
293  //Below from all maximum energy CaloParticles, we save the one with the largest amount
294  //of related recHits.
295  for (auto& c : occurrencesCPinMCL) {
296  if (c.second > maxCPNumberOfHitsInMCL) {
297  maxCPId_byNumberOfHits = c.first;
298  maxCPNumberOfHitsInMCL = c.second;
299  }
300  }
301 
302  //Find the CaloParticle that has the maximum energy shared with the multiCluster under study.
303  for (auto& c : CPEnergyInMCL) {
304  if (c.second > maxEnergySharedMCLandCP) {
305  maxCPId_byEnergy = c.first;
306  maxEnergySharedMCLandCP = c.second;
307  }
308  }
309  //The energy of the CaloParticle that found to have the maximum energy shared with the multiCluster under study.
310  float totalCPEnergyFromLayerCP = 0.f;
311  if (maxCPId_byEnergy >= 0) {
312  //Loop through all layers
313  for (unsigned int j = 0; j < layers_ * 2; ++j) {
314  totalCPEnergyFromLayerCP = totalCPEnergyFromLayerCP + cPOnLayer[maxCPId_byEnergy][j].energy;
315  }
316  energyFractionOfCPinMCL = maxEnergySharedMCLandCP / totalCPEnergyFromLayerCP;
317  if (mClusters[mcId].energy() > 0.f) {
318  energyFractionOfMCLinCP = maxEnergySharedMCLandCP / mClusters[mcId].energy();
319  }
320  }
321 
322  LogDebug("MultiClusterAssociatorByEnergyScoreImpl") << std::setw(12) << "multiCluster"
323  << "\t" << std::setw(10) << "mulcl energy"
324  << "\t" << std::setw(5) << "nhits"
325  << "\t" << std::setw(12) << "noise hits"
326  << "\t" << std::setw(22) << "maxCPId_byNumberOfHits"
327  << "\t" << std::setw(8) << "nhitsCP"
328  << "\t" << std::setw(16) << "maxCPId_byEnergy"
329  << "\t" << std::setw(23) << "maxEnergySharedMCLandCP"
330  << "\t" << std::setw(22) << "totalCPEnergyFromAllLayerCP"
331  << "\t" << std::setw(22) << "energyFractionOfMCLinCP"
332  << "\t" << std::setw(25) << "energyFractionOfCPinMCL"
333  << "\t" << std::endl;
334  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
335  << std::setw(12) << mcId << "\t" << std::setw(10) << mClusters[mcId].energy() << "\t" << std::setw(5)
336  << numberOfHitsInMCL << "\t" << std::setw(12) << numberOfNoiseHitsInMCL << "\t" << std::setw(22)
337  << maxCPId_byNumberOfHits << "\t" << std::setw(8) << maxCPNumberOfHitsInMCL << "\t" << std::setw(16)
338  << maxCPId_byEnergy << "\t" << std::setw(23) << maxEnergySharedMCLandCP << "\t" << std::setw(22)
339  << totalCPEnergyFromLayerCP << "\t" << std::setw(22) << energyFractionOfMCLinCP << "\t" << std::setw(25)
340  << energyFractionOfCPinMCL << std::endl;
341  }
342  } // end of loop through multiClusters
343 
344  // Update cpsInMultiCluster; compute the score MultiCluster-to-CaloParticle,
345  // together with the returned AssociationMap
346  for (unsigned int mcId = 0; mcId < nMultiClusters; ++mcId) {
347  // find the unique caloParticles id contributing to the multilusters
348  std::sort(cpsInMultiCluster[mcId].begin(), cpsInMultiCluster[mcId].end());
349  auto last = std::unique(cpsInMultiCluster[mcId].begin(), cpsInMultiCluster[mcId].end());
350  cpsInMultiCluster[mcId].erase(last, cpsInMultiCluster[mcId].end());
351 
352  const auto& hits_and_fractions = mClusters[mcId].hitsAndFractions();
353  unsigned int numberOfHitsInLC = hits_and_fractions.size();
354  if (numberOfHitsInLC > 0) {
355  if (mClusters[mcId].energy() == 0. && !cpsInMultiCluster[mcId].empty()) {
356  //Loop through all CaloParticles contributing to multiCluster mcId.
357  for (auto& cpPair : cpsInMultiCluster[mcId]) {
358  //In case of a multiCluster with zero energy but related CaloParticles the score is set to 1.
359  cpPair.second = 1.;
360  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
361  << "multiClusterId : \t " << mcId << "\t CP id : \t" << cpPair.first << "\t score \t " << cpPair.second
362  << "\n";
363  }
364  continue;
365  }
366 
367  // Compute the correct normalization
368  float invMultiClusterEnergyWeight = 0.f;
369  for (auto const& haf : mClusters[mcId].hitsAndFractions()) {
370  invMultiClusterEnergyWeight +=
371  (haf.second * hitMap_->at(haf.first)->energy()) * (haf.second * hitMap_->at(haf.first)->energy());
372  }
373  invMultiClusterEnergyWeight = 1.f / invMultiClusterEnergyWeight;
374 
375  for (unsigned int i = 0; i < numberOfHitsInLC; ++i) {
376  DetId rh_detid = hits_and_fractions[i].first;
377  float rhFraction = hits_and_fractions[i].second;
378 
379  bool hitWithNoCP = (detIdToCaloParticleId_Map.find(rh_detid) == detIdToCaloParticleId_Map.end());
380 
381  auto itcheck = hitMap_->find(rh_detid);
382  const HGCRecHit* hit = itcheck->second;
383  float hitEnergyWeight = hit->energy() * hit->energy();
384 
385  for (auto& cpPair : cpsInMultiCluster[mcId]) {
386  unsigned int multiClusterId = cpPair.first;
387  float cpFraction = 0.f;
388  if (!hitWithNoCP) {
389  auto findHitIt = std::find(detIdToCaloParticleId_Map[rh_detid].begin(),
390  detIdToCaloParticleId_Map[rh_detid].end(),
391  hgcal::detIdInfoInCluster{multiClusterId, 0.f});
392  if (findHitIt != detIdToCaloParticleId_Map[rh_detid].end())
393  cpFraction = findHitIt->fraction;
394  }
395  if (cpPair.second == FLT_MAX) {
396  cpPair.second = 0.f;
397  }
398  cpPair.second +=
399  (rhFraction - cpFraction) * (rhFraction - cpFraction) * hitEnergyWeight * invMultiClusterEnergyWeight;
400  }
401  } // End of loop over Hits within a MultiCluster
402 #ifdef EDM_ML_DEBUG
403  //In case of a multiCluster with some energy but none related CaloParticles print some info.
404  if (cpsInMultiCluster[mcId].empty())
405  LogDebug("MultiClusterAssociatorByEnergyScoreImpl") << "multiCluster Id: \t" << mcId << "\tCP id:\t-1 "
406  << "\t score \t-1"
407  << "\n";
408 #endif
409  }
410  } // End of loop over MultiClusters
411 
412  // Compute the CaloParticle-To-MultiCluster score
413  for (const auto& cpId : cPSelectedIndices) {
414  for (unsigned int layerId = 0; layerId < layers_ * 2; ++layerId) {
415  unsigned int CPNumberOfHits = cPOnLayer[cpId][layerId].hits_and_fractions.size();
416  if (CPNumberOfHits == 0)
417  continue;
418 #ifdef EDM_ML_DEBUG
419  int mcWithMaxEnergyInCP = -1;
420  float maxEnergyMCLperlayerinCP = 0.f;
421  float CPenergy = cPOnLayer[cpId][layerId].energy;
422  float CPEnergyFractionInMCLperlayer = 0.f;
423  for (auto& mc : cPOnLayer[cpId][layerId].multiClusterIdToEnergyAndScore) {
424  if (mc.second.first > maxEnergyMCLperlayerinCP) {
425  maxEnergyMCLperlayerinCP = mc.second.first;
426  mcWithMaxEnergyInCP = mc.first;
427  }
428  }
429  if (CPenergy > 0.f)
430  CPEnergyFractionInMCLperlayer = maxEnergyMCLperlayerinCP / CPenergy;
431 
432  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
433  << std::setw(8) << "LayerId:\t" << std::setw(12) << "caloparticle\t" << std::setw(15) << "cp total energy\t"
434  << std::setw(15) << "cpEnergyOnLayer\t" << std::setw(14) << "CPNhitsOnLayer\t" << std::setw(18)
435  << "mcWithMaxEnergyInCP\t" << std::setw(15) << "maxEnergyMCLinCP\t" << std::setw(20)
436  << "CPEnergyFractionInMCL"
437  << "\n";
438  LogDebug("MultiClusterAssociatorByEnergyScoreImpl")
439  << std::setw(8) << layerId << "\t" << std::setw(12) << cpId << "\t" << std::setw(15)
440  << caloParticles[cpId].energy() << "\t" << std::setw(15) << CPenergy << "\t" << std::setw(14)
441  << CPNumberOfHits << "\t" << std::setw(18) << mcWithMaxEnergyInCP << "\t" << std::setw(15)
442  << maxEnergyMCLperlayerinCP << "\t" << std::setw(20) << CPEnergyFractionInMCLperlayer << "\n";
443 #endif
444 
445  for (unsigned int i = 0; i < CPNumberOfHits; ++i) {
446  auto& cp_hitDetId = cPOnLayer[cpId][layerId].hits_and_fractions[i].first;
447  auto& cpFraction = cPOnLayer[cpId][layerId].hits_and_fractions[i].second;
448 
449  bool hitWithNoMCL = false;
450  if (cpFraction == 0.f)
451  continue; //hopefully this should never happen
452  auto hit_find_in_MCL = detIdToMultiClusterId_Map.find(cp_hitDetId);
453  if (hit_find_in_MCL == detIdToMultiClusterId_Map.end())
454  hitWithNoMCL = true;
455  auto itcheck = hitMap_->find(cp_hitDetId);
456  const HGCRecHit* hit = itcheck->second;
457  float hitEnergyWeight = hit->energy() * hit->energy();
458  for (auto& mcPair : cPOnLayer[cpId][layerId].multiClusterIdToEnergyAndScore) {
459  unsigned int multiClusterId = mcPair.first;
460  float mcFraction = 0.f;
461 
462  if (!hitWithNoMCL) {
463  auto findHitIt = std::find(detIdToMultiClusterId_Map[cp_hitDetId].begin(),
464  detIdToMultiClusterId_Map[cp_hitDetId].end(),
465  hgcal::detIdInfoInMultiCluster{multiClusterId, 0, 0.f});
466  if (findHitIt != detIdToMultiClusterId_Map[cp_hitDetId].end())
467  mcFraction = findHitIt->fraction;
468  }
469  //Observe here that we do not divide as before by the layer cluster energy weight. We should sum first
470  //over all layers and divide with the total CP energy over all layers.
471  if (mcPair.second.second == FLT_MAX) {
472  mcPair.second.second = 0.f;
473  }
474  mcPair.second.second += (mcFraction - cpFraction) * (mcFraction - cpFraction) * hitEnergyWeight;
475 #ifdef EDM_ML_DEBUG
476  LogDebug("HGCalValidator") << "multiClusterId:\t" << multiClusterId << "\tmcfraction,cpfraction:\t"
477  << mcFraction << ", " << cpFraction << "\thitEnergyWeight:\t" << hitEnergyWeight
478  << "\tcurrent score numerator:\t" << mcPair.second.second << "\n";
479 #endif
480  } // End of loop over MultiClusters linked to hits of this CaloParticle
481  } // End of loop over hits of CaloParticle on a Layer
482 #ifdef EDM_ML_DEBUG
483  if (cPOnLayer[cpId][layerId].multiClusterIdToEnergyAndScore.empty())
484  LogDebug("HGCalValidator") << "CP Id: \t" << cpId << "\t MCL id:\t-1 "
485  << "\t layer \t " << layerId << " Sub score in \t -1"
486  << "\n";
487 
488 #endif
489  }
490  }
491  return {cpsInMultiCluster, cPOnLayer};
492 }
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
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:184
std::vector< std::vector< hgcal::caloParticleOnALayer > > caloParticleToMultiCluster
Monte Carlo truth information used for tracking validation.
Definition: SimCluster.h:29
def unique(seq, keepstr=True)
Definition: tier0.py:24
double f[11][100]
Definition: DetId.h:17
static void removeCPFromPU(const std::vector< CaloParticle > &caloParticles, std::vector< size_t > &cPIndices, bool hardScatterOnly=true)
const std::unordered_map< DetId, const HGCRecHit * > * hitMap_
std::vector< std::vector< std::pair< unsigned int, float > > > multiClusterToCaloParticle
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
def cp(fromDir, toDir, listOfFiles, overwrite=False, smallList=False)
#define LogDebug(id)

Member Data Documentation

◆ hardScatterOnly_

const bool MultiClusterAssociatorByEnergyScoreImpl::hardScatterOnly_
private

Definition at line 63 of file MultiClusterAssociatorByEnergyScoreImpl.h.

◆ hitMap_

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

Definition at line 65 of file MultiClusterAssociatorByEnergyScoreImpl.h.

Referenced by makeConnections().

◆ layers_

unsigned MultiClusterAssociatorByEnergyScoreImpl::layers_
private

◆ productGetter_

edm::EDProductGetter const* MultiClusterAssociatorByEnergyScoreImpl::productGetter_
private

◆ recHitTools_

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