Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-16 05:06:39

0001 // Author: Felice Pantaleo, felice.pantaleo@cern.ch 06/2024
0002 #include "TracksterToSimTracksterAssociatorByHitsProducer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "FWCore/Utilities/interface/EDGetToken.h"
0009 #include "DataFormats/HGCalReco/interface/Trackster.h"
0010 #include "SimDataFormats/Associations/interface/TICLAssociationMap.h"
0011 #include "DataFormats/Provenance/interface/ProductID.h"
0012 #include "DataFormats/HGCRecHit/interface/HGCRecHitCollections.h"
0013 #include "CommonTools/RecoAlgos/interface/MultiVectorManager.h"
0014 #include "SimDataFormats/CaloAnalysis/interface/CaloParticle.h"
0015 #include "SimDataFormats/CaloAnalysis/interface/SimCluster.h"
0016 
0017 TracksterToSimTracksterAssociatorByHitsProducer::TracksterToSimTracksterAssociatorByHitsProducer(
0018     const edm::ParameterSet& pset)
0019     : recoTracksterCollectionToken_(
0020           consumes<std::vector<ticl::Trackster>>(pset.getParameter<edm::InputTag>("tracksters"))),
0021       simTracksterCollectionToken_(
0022           consumes<std::vector<ticl::Trackster>>(pset.getParameter<edm::InputTag>("simTracksters"))),
0023       simTracksterFromCPCollectionToken_(
0024           consumes<std::vector<ticl::Trackster>>(pset.getParameter<edm::InputTag>("simTrackstersFromCP"))),
0025       hitToTracksterMapToken_(
0026           consumes<ticl::AssociationMap<ticl::mapWithFraction>>(pset.getParameter<edm::InputTag>("hitToTracksterMap"))),
0027       hitToSimTracksterMapToken_(consumes<ticl::AssociationMap<ticl::mapWithFraction>>(
0028           pset.getParameter<edm::InputTag>("hitToSimTracksterMap"))),
0029       hitToSimTracksterFromCPMapToken_(consumes<ticl::AssociationMap<ticl::mapWithFraction>>(
0030           pset.getParameter<edm::InputTag>("hitToSimTracksterFromCPMap"))),
0031       hitToSimClusterMapToken_(consumes<ticl::AssociationMap<ticl::mapWithFraction>>(
0032           pset.getParameter<edm::InputTag>("hitToSimClusterMap"))),
0033       hitToCaloParticleMapToken_(consumes<ticl::AssociationMap<ticl::mapWithFraction>>(
0034           pset.getParameter<edm::InputTag>("hitToCaloParticleMap"))),
0035       tracksterToHitMapToken_(
0036           consumes<ticl::AssociationMap<ticl::mapWithFraction>>(pset.getParameter<edm::InputTag>("tracksterToHitMap"))),
0037       simTracksterToHitMapToken_(consumes<ticl::AssociationMap<ticl::mapWithFraction>>(
0038           pset.getParameter<edm::InputTag>("simTracksterToHitMap"))),
0039       simTracksterFromCPToHitMapToken_(consumes<ticl::AssociationMap<ticl::mapWithFraction>>(
0040           pset.getParameter<edm::InputTag>("simTracksterFromCPToHitMap"))),
0041       caloParticleToken_(consumes<std::vector<CaloParticle>>(pset.getParameter<edm::InputTag>("caloParticles"))) {
0042   auto hitsTags = pset.getParameter<std::vector<edm::InputTag>>("hits");
0043   for (const auto& tag : hitsTags) {
0044     hitsTokens_.push_back(consumes<HGCRecHitCollection>(tag));
0045   }
0046   produces<
0047       ticl::AssociationMap<ticl::mapWithFractionAndScore, std::vector<ticl::Trackster>, std::vector<ticl::Trackster>>>(
0048       "tracksterToSimTracksterMap");
0049   produces<
0050       ticl::AssociationMap<ticl::mapWithFractionAndScore, std::vector<ticl::Trackster>, std::vector<ticl::Trackster>>>(
0051       "simTracksterToTracksterMap");
0052   produces<
0053       ticl::AssociationMap<ticl::mapWithFractionAndScore, std::vector<ticl::Trackster>, std::vector<ticl::Trackster>>>(
0054       "tracksterToSimTracksterFromCPMap");
0055   produces<
0056       ticl::AssociationMap<ticl::mapWithFractionAndScore, std::vector<ticl::Trackster>, std::vector<ticl::Trackster>>>(
0057       "simTracksterFromCPToTracksterMap");
0058 }
0059 
0060 TracksterToSimTracksterAssociatorByHitsProducer::~TracksterToSimTracksterAssociatorByHitsProducer() {}
0061 
0062 void TracksterToSimTracksterAssociatorByHitsProducer::produce(edm::StreamID,
0063                                                               edm::Event& iEvent,
0064                                                               const edm::EventSetup& iSetup) const {
0065   using namespace edm;
0066 
0067   Handle<std::vector<ticl::Trackster>> recoTrackstersHandle;
0068   iEvent.getByToken(recoTracksterCollectionToken_, recoTrackstersHandle);
0069   const auto& recoTracksters = *recoTrackstersHandle;
0070 
0071   Handle<std::vector<ticl::Trackster>> simTrackstersHandle;
0072   iEvent.getByToken(simTracksterCollectionToken_, simTrackstersHandle);
0073   const auto& simTracksters = *simTrackstersHandle;
0074 
0075   Handle<std::vector<ticl::Trackster>> simTrackstersFromCPHandle;
0076   iEvent.getByToken(simTracksterFromCPCollectionToken_, simTrackstersFromCPHandle);
0077   const auto& simTrackstersFromCP = *simTrackstersFromCPHandle;
0078 
0079   Handle<ticl::AssociationMap<ticl::mapWithFraction>> hitToTracksterMapHandle;
0080   iEvent.getByToken(hitToTracksterMapToken_, hitToTracksterMapHandle);
0081   const auto& hitToTracksterMap = *hitToTracksterMapHandle;
0082 
0083   Handle<ticl::AssociationMap<ticl::mapWithFraction>> hitToSimTracksterMapHandle;
0084   iEvent.getByToken(hitToSimTracksterMapToken_, hitToSimTracksterMapHandle);
0085   const auto& hitToSimTracksterMap = *hitToSimTracksterMapHandle;
0086 
0087   Handle<ticl::AssociationMap<ticl::mapWithFraction>> hitToSimTracksterFromCPMapHandle;
0088   iEvent.getByToken(hitToSimTracksterFromCPMapToken_, hitToSimTracksterFromCPMapHandle);
0089   const auto& hitToSimTracksterFromCPMap = *hitToSimTracksterFromCPMapHandle;
0090 
0091   Handle<ticl::AssociationMap<ticl::mapWithFraction>> hitToSimClusterMapHandle;
0092   iEvent.getByToken(hitToSimClusterMapToken_, hitToSimClusterMapHandle);
0093   const auto& hitToSimClusterMap = *hitToSimClusterMapHandle;
0094 
0095   Handle<ticl::AssociationMap<ticl::mapWithFraction>> hitToCaloParticleMapHandle;
0096   iEvent.getByToken(hitToCaloParticleMapToken_, hitToCaloParticleMapHandle);
0097   const auto& hitToCaloParticleMap = *hitToCaloParticleMapHandle;
0098 
0099   Handle<ticl::AssociationMap<ticl::mapWithFraction>> tracksterToHitMapHandle;
0100   iEvent.getByToken(tracksterToHitMapToken_, tracksterToHitMapHandle);
0101   const auto& tracksterToHitMap = *tracksterToHitMapHandle;
0102 
0103   Handle<ticl::AssociationMap<ticl::mapWithFraction>> simTracksterToHitMapHandle;
0104   iEvent.getByToken(simTracksterToHitMapToken_, simTracksterToHitMapHandle);
0105   const auto& simTracksterToHitMap = *simTracksterToHitMapHandle;
0106 
0107   Handle<ticl::AssociationMap<ticl::mapWithFraction>> simTracksterFromCPToHitMapHandle;
0108   iEvent.getByToken(simTracksterFromCPToHitMapToken_, simTracksterFromCPToHitMapHandle);
0109   const auto& simTracksterFromCPToHitMap = *simTracksterFromCPToHitMapHandle;
0110 
0111   Handle<std::vector<CaloParticle>> caloParticlesHandle;
0112   iEvent.getByToken(caloParticleToken_, caloParticlesHandle);
0113 
0114   MultiVectorManager<HGCRecHit> rechitManager;
0115   for (const auto& token : hitsTokens_) {
0116     Handle<HGCRecHitCollection> hitsHandle;
0117     iEvent.getByToken(token, hitsHandle);
0118     rechitManager.addVector(*hitsHandle);
0119   }
0120 
0121   auto tracksterToSimTracksterMap = std::make_unique<
0122       ticl::AssociationMap<ticl::mapWithFractionAndScore, std::vector<ticl::Trackster>, std::vector<ticl::Trackster>>>(
0123       recoTrackstersHandle, simTrackstersHandle, iEvent);
0124   auto tracksterToSimTracksterFromCPMap = std::make_unique<
0125       ticl::AssociationMap<ticl::mapWithFractionAndScore, std::vector<ticl::Trackster>, std::vector<ticl::Trackster>>>(
0126       recoTrackstersHandle, simTrackstersFromCPHandle, iEvent);
0127 
0128   auto simTracksterToTracksterMap = std::make_unique<
0129       ticl::AssociationMap<ticl::mapWithFractionAndScore, std::vector<ticl::Trackster>, std::vector<ticl::Trackster>>>(
0130       simTrackstersHandle, recoTrackstersHandle, iEvent);
0131   auto simTracksterFromCPToTracksterMap = std::make_unique<
0132       ticl::AssociationMap<ticl::mapWithFractionAndScore, std::vector<ticl::Trackster>, std::vector<ticl::Trackster>>>(
0133       simTrackstersFromCPHandle, recoTrackstersHandle, iEvent);
0134   for (unsigned int tracksterIndex = 0; tracksterIndex < recoTracksters.size(); ++tracksterIndex) {
0135     edm::Ref<std::vector<ticl::Trackster>> recoTracksterRef(recoTrackstersHandle, tracksterIndex);
0136     float recoToSimScoresDenominator = 0.f;
0137     const auto& recoTracksterHitsAndFractions = tracksterToHitMap[tracksterIndex];
0138     ticl::AssociationMap<ticl::mapWithFraction> hitToAssociatedSimTracksterMap(recoTracksterHitsAndFractions.size());
0139     std::vector<unsigned int> associatedSimTracksterIndices;
0140     ticl::AssociationMap<ticl::mapWithFraction> hitToAssociatedSimTracksterFromCPMap(
0141         recoTracksterHitsAndFractions.size());
0142     std::vector<unsigned int> associatedSimTracksterFromCPIndices;
0143     for (unsigned int i = 0; i < recoTracksterHitsAndFractions.size(); ++i) {
0144       const auto& [hitIndex, recoFraction] = recoTracksterHitsAndFractions[i];
0145       const auto& recHit = rechitManager[hitIndex];
0146       float squaredRecoFraction = recoFraction * recoFraction;
0147       float rechitEnergy = recHit.energy();
0148       float squaredRecHitEnergy = rechitEnergy * rechitEnergy;
0149       recoToSimScoresDenominator += squaredRecoFraction * squaredRecHitEnergy;
0150 
0151       const auto& hitToSimTracksterVec = hitToSimTracksterMap[hitIndex];
0152       for (const auto& [simTracksterIndex, fraction] : hitToSimTracksterVec) {
0153         const auto& simTrackster = simTracksters[simTracksterIndex];
0154         auto& seed = simTrackster.seedID();
0155         float simFraction = 0;
0156         if (seed == caloParticlesHandle.id()) {
0157           unsigned int caloParticleIndex = simTrackster.seedIndex();
0158           auto it = std::find_if(hitToCaloParticleMap[hitIndex].begin(),
0159                                  hitToCaloParticleMap[hitIndex].end(),
0160                                  [caloParticleIndex](const auto& pair) { return pair.first == caloParticleIndex; });
0161           if (it != hitToCaloParticleMap[hitIndex].end()) {
0162             simFraction = it->second;
0163           }
0164         } else {
0165           unsigned int simClusterIndex = simTracksters[simTracksterIndex].seedIndex();
0166           auto it = std::find_if(hitToSimClusterMap[hitIndex].begin(),
0167                                  hitToSimClusterMap[hitIndex].end(),
0168                                  [simClusterIndex](const auto& pair) { return pair.first == simClusterIndex; });
0169           if (it != hitToSimClusterMap[hitIndex].end()) {
0170             simFraction = it->second;
0171           }
0172         }
0173         hitToAssociatedSimTracksterMap.insert(i, simTracksterIndex, simFraction);
0174         associatedSimTracksterIndices.push_back(simTracksterIndex);
0175       }
0176 
0177       // do the same for caloparticles and simTracksterFromCP
0178       const auto& hitToSimTracksterFromCPVec = hitToSimTracksterFromCPMap[hitIndex];
0179       for (const auto& [simTracksterIndex, simFraction] : hitToSimTracksterFromCPVec) {
0180         unsigned int caloParticleIndex = simTracksters[simTracksterIndex].seedIndex();
0181         float caloParticleFraction = 0;
0182         auto it = std::find_if(hitToCaloParticleMap[hitIndex].begin(),
0183                                hitToCaloParticleMap[hitIndex].end(),
0184                                [caloParticleIndex](const auto& pair) { return pair.first == caloParticleIndex; });
0185         if (it != hitToCaloParticleMap[hitIndex].end()) {
0186           caloParticleFraction = it->second;
0187         }
0188         hitToAssociatedSimTracksterFromCPMap.insert(i, simTracksterIndex, caloParticleFraction);
0189         associatedSimTracksterFromCPIndices.push_back(simTracksterIndex);
0190       }
0191     }
0192     std::sort(associatedSimTracksterIndices.begin(), associatedSimTracksterIndices.end());
0193     associatedSimTracksterIndices.erase(
0194         std::unique(associatedSimTracksterIndices.begin(), associatedSimTracksterIndices.end()),
0195         associatedSimTracksterIndices.end());
0196 
0197     std::sort(associatedSimTracksterFromCPIndices.begin(), associatedSimTracksterFromCPIndices.end());
0198     associatedSimTracksterFromCPIndices.erase(
0199         std::unique(associatedSimTracksterFromCPIndices.begin(), associatedSimTracksterFromCPIndices.end()),
0200         associatedSimTracksterFromCPIndices.end());
0201 
0202     // Add missing sim tracksters with 0 shared energy to hitToAssociatedSimTracksterMap and hitToAssociatedSimTracksterFromCPMap
0203     for (unsigned int i = 0; i < recoTracksterHitsAndFractions.size(); ++i) {
0204       unsigned int hitId = recoTracksterHitsAndFractions[i].first;
0205       const auto& simTracksterVec = hitToSimTracksterMap[hitId];
0206       for (unsigned int simTracksterIndex : associatedSimTracksterIndices) {
0207         if (std::find_if(simTracksterVec.begin(), simTracksterVec.end(), [simTracksterIndex](const auto& pair) {
0208               return pair.first == simTracksterIndex;
0209             }) == simTracksterVec.end()) {
0210           hitToAssociatedSimTracksterMap.insert(i, simTracksterIndex, 0);
0211         }
0212       }
0213 
0214       const auto& simTracksterFromCPVec = hitToSimTracksterFromCPMap[hitId];
0215       for (unsigned int simTracksterIndex : associatedSimTracksterFromCPIndices) {
0216         if (std::find_if(
0217                 simTracksterFromCPVec.begin(), simTracksterFromCPVec.end(), [simTracksterIndex](const auto& pair) {
0218                   return pair.first == simTracksterIndex;
0219                 }) == simTracksterFromCPVec.end()) {
0220           hitToAssociatedSimTracksterFromCPMap.insert(i, simTracksterIndex, 0);
0221         }
0222       }
0223     }
0224 
0225     const float invDenominator = 1.f / recoToSimScoresDenominator;
0226 
0227     for (unsigned int i = 0; i < recoTracksterHitsAndFractions.size(); ++i) {
0228       unsigned int hitIndex = recoTracksterHitsAndFractions[i].first;
0229       const auto& recHit = rechitManager[hitIndex];
0230       float recoFraction = recoTracksterHitsAndFractions[i].second;
0231       float squaredRecoFraction = recoFraction * recoFraction;
0232       float squaredRecHitEnergy = recHit.energy() * recHit.energy();
0233       float recoSharedEnergy = recHit.energy() * recoFraction;
0234       const auto& simTracksterVec = hitToAssociatedSimTracksterMap[i];
0235       for (const auto& [simTracksterIndex, simFraction] : simTracksterVec) {
0236         edm::Ref<std::vector<ticl::Trackster>> simTracksterRef(simTrackstersHandle, simTracksterIndex);
0237         float sharedEnergy = std::min(simFraction * recHit.energy(), recoSharedEnergy);
0238         float squaredFraction =
0239             std::min(squaredRecoFraction, (recoFraction - simFraction) * (recoFraction - simFraction));
0240         float score = invDenominator * squaredFraction * squaredRecHitEnergy;
0241         tracksterToSimTracksterMap->insert(recoTracksterRef, simTracksterRef, sharedEnergy, score);
0242       }
0243 
0244       const auto& simTracksterFromCPVec = hitToAssociatedSimTracksterFromCPMap[i];
0245       for (const auto& [simTracksterIndex, simFraction] : simTracksterFromCPVec) {
0246         edm::Ref<std::vector<ticl::Trackster>> simTracksterRef(simTrackstersFromCPHandle, simTracksterIndex);
0247         float sharedEnergy = std::min(simFraction * recHit.energy(), recoSharedEnergy);
0248         float squaredFraction =
0249             std::min(squaredRecoFraction, (recoFraction - simFraction) * (recoFraction - simFraction));
0250         float score = invDenominator * squaredFraction * squaredRecHitEnergy;
0251         tracksterToSimTracksterFromCPMap->insert(recoTracksterRef, simTracksterRef, sharedEnergy, score);
0252       }
0253     }
0254   }
0255 
0256   // Reverse mapping: SimTrackster -> RecoTrackster
0257   for (unsigned int tracksterIndex = 0; tracksterIndex < simTracksters.size(); ++tracksterIndex) {
0258     edm::Ref<std::vector<ticl::Trackster>> simTracksterRef(simTrackstersHandle, tracksterIndex);
0259     float simToRecoScoresDenominator = 0.f;
0260     const auto& simTracksterHitsAndFractions = simTracksterToHitMap[tracksterIndex];
0261     ticl::AssociationMap<ticl::mapWithFraction> hitToAssociatedRecoTracksterMap(simTracksterHitsAndFractions.size());
0262     std::vector<unsigned int> associatedRecoTracksterIndices;
0263     const auto& simTrackster = simTracksters[tracksterIndex];
0264     auto& seed = simTrackster.seedID();
0265     unsigned int simObjectIndex = simTrackster.seedIndex();
0266     bool isSimTracksterFromCP = (seed == caloParticlesHandle.id());
0267     std::vector<float> simFractions(simTracksterHitsAndFractions.size(), 0.f);
0268     for (unsigned int i = 0; i < simTracksterHitsAndFractions.size(); ++i) {
0269       const auto& [hitIndex, simTracksterFraction] = simTracksterHitsAndFractions[i];
0270 
0271       auto it = isSimTracksterFromCP
0272                     ? (std::find_if(hitToCaloParticleMap[hitIndex].begin(),
0273                                     hitToCaloParticleMap[hitIndex].end(),
0274                                     [simObjectIndex](const auto& pair) { return pair.first == simObjectIndex; }))
0275                     : std::find_if(hitToSimClusterMap[hitIndex].begin(),
0276                                    hitToSimClusterMap[hitIndex].end(),
0277                                    [simObjectIndex](const auto& pair) { return pair.first == simObjectIndex; });
0278       if (it != hitToCaloParticleMap[hitIndex].end() and it != hitToSimClusterMap[hitIndex].end()) {
0279         simFractions[i] = it->second;
0280       }
0281       float simFraction = simFractions[i];
0282       const auto& recHit = rechitManager[hitIndex];
0283       float squaredSimFraction = simFraction * simFraction;
0284       float squaredRecHitEnergy = recHit.energy() * recHit.energy();
0285       simToRecoScoresDenominator += squaredSimFraction * squaredRecHitEnergy;
0286 
0287       const auto& hitToRecoTracksterVec = hitToTracksterMap[hitIndex];
0288       for (const auto& [recoTracksterIndex, recoFraction] : hitToRecoTracksterVec) {
0289         hitToAssociatedRecoTracksterMap.insert(i, recoTracksterIndex, recoFraction);
0290         associatedRecoTracksterIndices.push_back(recoTracksterIndex);
0291       }
0292     }
0293 
0294     std::sort(associatedRecoTracksterIndices.begin(), associatedRecoTracksterIndices.end());
0295     associatedRecoTracksterIndices.erase(
0296         std::unique(associatedRecoTracksterIndices.begin(), associatedRecoTracksterIndices.end()),
0297         associatedRecoTracksterIndices.end());
0298 
0299     for (unsigned int i = 0; i < simTracksterHitsAndFractions.size(); ++i) {
0300       unsigned int hitIndex = simTracksterHitsAndFractions[i].first;
0301       const auto& hitToRecoTracksterVec = hitToTracksterMap[hitIndex];
0302       for (unsigned int recoTracksterIndex : associatedRecoTracksterIndices) {
0303         if (std::find_if(
0304                 hitToRecoTracksterVec.begin(), hitToRecoTracksterVec.end(), [recoTracksterIndex](const auto& pair) {
0305                   return pair.first == recoTracksterIndex;
0306                 }) == hitToRecoTracksterVec.end()) {
0307           hitToAssociatedRecoTracksterMap.insert(i, recoTracksterIndex, 0);
0308         }
0309       }
0310     }
0311 
0312     const float invDenominator = 1.f / simToRecoScoresDenominator;
0313 
0314     for (unsigned int i = 0; i < simTracksterHitsAndFractions.size(); ++i) {
0315       const auto& [hitIndex, simTracksterFraction] = simTracksterHitsAndFractions[i];
0316       float simFraction = simFractions[i];
0317       const auto& recHit = rechitManager[hitIndex];
0318       float squaredSimFraction = simFraction * simFraction;
0319       float squaredRecHitEnergy = recHit.energy() * recHit.energy();
0320       float simSharedEnergy = recHit.energy() * simFraction;
0321 
0322       const auto& hitToRecoTracksterVec = hitToAssociatedRecoTracksterMap[i];
0323       for (const auto& [recoTracksterIndex, recoFraction] : hitToRecoTracksterVec) {
0324         edm::Ref<std::vector<ticl::Trackster>> recoTracksterRef(recoTrackstersHandle, recoTracksterIndex);
0325         float sharedEnergy = std::min(recoFraction * recHit.energy(), simSharedEnergy);
0326         float squaredFraction =
0327             std::min(squaredSimFraction, (recoFraction - simFraction) * (recoFraction - simFraction));
0328         float score = invDenominator * squaredFraction * squaredRecHitEnergy;
0329         simTracksterToTracksterMap->insert(simTracksterRef, recoTracksterRef, sharedEnergy, score);
0330       }
0331     }
0332   }
0333 
0334   // Repeat the reverse mapping process for SimTracksterFromCP
0335   for (unsigned int tracksterIndex = 0; tracksterIndex < simTrackstersFromCP.size(); ++tracksterIndex) {
0336     edm::Ref<std::vector<ticl::Trackster>> simTracksterRef(simTrackstersFromCPHandle, tracksterIndex);
0337     float simToRecoScoresDenominator = 0.f;
0338     const auto& simTracksterHitsAndFractions = simTracksterFromCPToHitMap[tracksterIndex];
0339     ticl::AssociationMap<ticl::mapWithFraction> hitToAssociatedRecoTracksterMap(simTracksterHitsAndFractions.size());
0340     std::vector<unsigned int> associatedRecoTracksterIndices;
0341     std::vector<float> simFractions(simTracksterHitsAndFractions.size(), 0.f);
0342     const auto& simTrackster = simTrackstersFromCP[tracksterIndex];
0343     unsigned int simObjectIndex = simTrackster.seedIndex();
0344     for (unsigned int i = 0; i < simTracksterHitsAndFractions.size(); ++i) {
0345       const auto& [hitIndex, simTracksterFraction] = simTracksterHitsAndFractions[i];
0346       auto it = std::find_if(hitToCaloParticleMap[hitIndex].begin(),
0347                              hitToCaloParticleMap[hitIndex].end(),
0348                              [simObjectIndex](const auto& pair) { return pair.first == simObjectIndex; });
0349       if (it != hitToCaloParticleMap[hitIndex].end()) {
0350         simFractions[i] = it->second;
0351       }
0352       float simFraction = simFractions[i];
0353 
0354       const auto& recHit = rechitManager[hitIndex];
0355       float squaredSimFraction = simFraction * simFraction;
0356       float squaredRecHitEnergy = recHit.energy() * recHit.energy();
0357       simToRecoScoresDenominator += squaredSimFraction * squaredRecHitEnergy;
0358 
0359       const auto& hitToRecoTracksterVec = hitToTracksterMap[hitIndex];
0360       for (const auto& [recoTracksterIndex, recoFraction] : hitToRecoTracksterVec) {
0361         hitToAssociatedRecoTracksterMap.insert(i, recoTracksterIndex, recoFraction);
0362         associatedRecoTracksterIndices.push_back(recoTracksterIndex);
0363       }
0364     }
0365 
0366     std::sort(associatedRecoTracksterIndices.begin(), associatedRecoTracksterIndices.end());
0367     associatedRecoTracksterIndices.erase(
0368         std::unique(associatedRecoTracksterIndices.begin(), associatedRecoTracksterIndices.end()),
0369         associatedRecoTracksterIndices.end());
0370 
0371     for (unsigned int i = 0; i < simTracksterHitsAndFractions.size(); ++i) {
0372       unsigned int hitIndex = simTracksterHitsAndFractions[i].first;
0373       const auto& hitToRecoTracksterVec = hitToTracksterMap[hitIndex];
0374       for (unsigned int recoTracksterIndex : associatedRecoTracksterIndices) {
0375         if (std::find_if(
0376                 hitToRecoTracksterVec.begin(), hitToRecoTracksterVec.end(), [recoTracksterIndex](const auto& pair) {
0377                   return pair.first == recoTracksterIndex;
0378                 }) == hitToRecoTracksterVec.end()) {
0379           hitToAssociatedRecoTracksterMap.insert(i, recoTracksterIndex, 0.f);
0380         }
0381       }
0382     }
0383 
0384     const float invDenominator = 1.f / simToRecoScoresDenominator;
0385 
0386     for (unsigned int i = 0; i < simTracksterHitsAndFractions.size(); ++i) {
0387       const auto& [hitIndex, simTracksterFraction] = simTracksterHitsAndFractions[i];
0388       const auto& recHit = rechitManager[hitIndex];
0389       float simFraction = simFractions[i];
0390       float squaredSimFraction = simFraction * simFraction;
0391       float squaredRecHitEnergy = recHit.energy() * recHit.energy();
0392       float simSharedEnergy = recHit.energy() * simFraction;
0393 
0394       const auto& hitToRecoTracksterVec = hitToAssociatedRecoTracksterMap[i];
0395       for (const auto& [recoTracksterIndex, recoFraction] : hitToRecoTracksterVec) {
0396         edm::Ref<std::vector<ticl::Trackster>> recoTracksterRef(recoTrackstersHandle, recoTracksterIndex);
0397         float sharedEnergy = std::min(recoFraction * recHit.energy(), simSharedEnergy);
0398         float squaredFraction =
0399             std::min(squaredSimFraction, (recoFraction - simFraction) * (recoFraction - simFraction));
0400         float score = invDenominator * squaredFraction * squaredRecHitEnergy;
0401         simTracksterFromCPToTracksterMap->insert(simTracksterRef, recoTracksterRef, sharedEnergy, score);
0402       }
0403     }
0404   }
0405   tracksterToSimTracksterMap->sort(true);
0406   tracksterToSimTracksterFromCPMap->sort(true);
0407   simTracksterToTracksterMap->sort(true);
0408   simTracksterFromCPToTracksterMap->sort(true);
0409 
0410   iEvent.put(std::move(tracksterToSimTracksterMap), "tracksterToSimTracksterMap");
0411   iEvent.put(std::move(tracksterToSimTracksterFromCPMap), "tracksterToSimTracksterFromCPMap");
0412   iEvent.put(std::move(simTracksterToTracksterMap), "simTracksterToTracksterMap");
0413   iEvent.put(std::move(simTracksterFromCPToTracksterMap), "simTracksterFromCPToTracksterMap");
0414 }
0415 
0416 void TracksterToSimTracksterAssociatorByHitsProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0417   edm::ParameterSetDescription desc;
0418   desc.add<edm::InputTag>("tracksters", edm::InputTag("ticlTrackstersMerge"));
0419   desc.add<edm::InputTag>("simTracksters", edm::InputTag("ticlSimTracksters"));
0420   desc.add<edm::InputTag>("simTrackstersFromCP", edm::InputTag("ticlSimTracksters", "fromCPs"));
0421 
0422   desc.add<edm::InputTag>("hitToTracksterMap", edm::InputTag("hitToTracksterAssociator", "hitToTracksterMap"));
0423   desc.add<edm::InputTag>("hitToSimTracksterMap",
0424                           edm::InputTag("allHitToTracksterAssociations", "hitToticlSimTracksters"));
0425   desc.add<edm::InputTag>("hitToSimTracksterFromCPMap",
0426                           edm::InputTag("allHitToTracksterAssociations", "hitToticlSimTrackstersfromCPs"));
0427   desc.add<edm::InputTag>("hitToSimClusterMap",
0428                           edm::InputTag("hitToSimClusterCaloParticleAssociator", "hitToSimClusterMap"));
0429   desc.add<edm::InputTag>("hitToCaloParticleMap",
0430                           edm::InputTag("hitToSimClusterCaloParticleAssociator", "hitToCaloParticleMap"));
0431   desc.add<edm::InputTag>("tracksterToHitMap", edm::InputTag("hitToTrackstersAssociationPR", "tracksterToHitMap"));
0432   desc.add<edm::InputTag>("simTracksterToHitMap",
0433                           edm::InputTag("allHitToTracksterAssociations", "ticlSimTrackstersToHit"));
0434   desc.add<edm::InputTag>("simTracksterFromCPToHitMap",
0435                           edm::InputTag("allHitToTracksterAssociations", "ticlSimTrackstersfromCPsToHit"));
0436   desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
0437 
0438   desc.add<std::vector<edm::InputTag>>("hits",
0439                                        {edm::InputTag("HGCalRecHit", "HGCEERecHits"),
0440                                         edm::InputTag("HGCalRecHit", "HGCHEFRecHits"),
0441                                         edm::InputTag("HGCalRecHit", "HGCHEBRecHits")});
0442   descriptions.add("tracksterToSimTracksterAssociatorByHitsProducer", desc);
0443 }
0444 
0445 // Define this as a plug-in
0446 DEFINE_FWK_MODULE(TracksterToSimTracksterAssociatorByHitsProducer);