File indexing completed on 2024-09-07 04:38:11
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "SimTracker/TrackTriggerAssociation/plugins/TTStubAssociator.h"
0011
0012
0013 template <>
0014 void TTStubAssociator<Ref_Phase2TrackerDigi_>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0015
0016 if (iEvent.isRealData())
0017 return;
0018
0019
0020 if (ttClusterTruthInputTags_.size() != ttStubsInputTags_.size()) {
0021 edm::LogError("TTStubAsso ") << "E R R O R! the InputTag vectors have different size!";
0022 return;
0023 }
0024
0025 int ncont1 = 0;
0026
0027 const TrackerGeometry* const theTrackerGeom = &iSetup.getData(theTrackerGeometryToken_);
0028 const TrackerTopology* const tTopo = &iSetup.getData(theTrackerTopologyToken_);
0029
0030
0031
0032 for (const auto& iTag : ttStubsTokens_) {
0033
0034 auto associationMapForOutput = std::make_unique<TTStubAssociationMap<Ref_Phase2TrackerDigi_>>();
0035
0036
0037 edm::Handle<TTStubDetSetVec> ttStubHandle;
0038 iEvent.getByToken(iTag, ttStubHandle);
0039
0040
0041 edm::Handle<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> ttClusterAssociationMapHandle;
0042 iEvent.getByToken(ttClusterTruthTokens_.at(ncont1), ttClusterAssociationMapHandle);
0043
0044
0045 std::map<TTStubRef, TrackingParticlePtr> stubToTrackingParticleMap;
0046 std::map<TrackingParticlePtr, std::vector<TTStubRef>> trackingParticleToStubVectorMap;
0047
0048
0049
0050 if (not ttStubHandle->empty()) {
0051 for (const auto& gd : theTrackerGeom->dets()) {
0052 DetId detid = gd->geographicalId();
0053 if (detid.subdetId() != StripSubdetector::TOB && detid.subdetId() != StripSubdetector::TID)
0054 continue;
0055
0056 if (!tTopo->isLower(detid))
0057 continue;
0058
0059 DetId stackDetid = tTopo->stack(detid);
0060
0061 if (ttStubHandle->find(stackDetid) == ttStubHandle->end())
0062 continue;
0063
0064
0065 edmNew::DetSet<TTStub<Ref_Phase2TrackerDigi_>> stubs = (*ttStubHandle)[stackDetid];
0066
0067 for (auto contentIter = stubs.begin(); contentIter != stubs.end(); ++contentIter) {
0068
0069 TTStubRef tempStubRef = edmNew::makeRefTo(ttStubHandle, contentIter);
0070
0071
0072 for (unsigned int ic = 0; ic < 2; ic++) {
0073 const std::vector<TrackingParticlePtr>& tempTPs =
0074 ttClusterAssociationMapHandle->findTrackingParticlePtrs(tempStubRef->clusterRef(ic));
0075
0076 for (const TrackingParticlePtr& testTP : tempTPs) {
0077 if (testTP.isNull())
0078 continue;
0079
0080
0081 if (trackingParticleToStubVectorMap.find(testTP) == trackingParticleToStubVectorMap.end()) {
0082 std::vector<TTStubRef> stubVector;
0083 trackingParticleToStubVectorMap.emplace(testTP, stubVector);
0084 }
0085 trackingParticleToStubVectorMap.find(testTP)->second.push_back(tempStubRef);
0086 }
0087 }
0088
0089
0090
0091
0092 if (ttClusterAssociationMapHandle->isUnknown(tempStubRef->clusterRef(0)) ||
0093 ttClusterAssociationMapHandle->isUnknown(tempStubRef->clusterRef(1))) {
0094
0095
0096
0097 continue;
0098 } else {
0099
0100
0101
0102 if (ttClusterAssociationMapHandle->isGenuine(tempStubRef->clusterRef(0)) &&
0103 ttClusterAssociationMapHandle->isGenuine(tempStubRef->clusterRef(1))) {
0104
0105
0106
0107 if (ttClusterAssociationMapHandle->findTrackingParticlePtr(tempStubRef->clusterRef(0)).get() ==
0108 ttClusterAssociationMapHandle->findTrackingParticlePtr(tempStubRef->clusterRef(1)).get()) {
0109
0110 const TrackingParticlePtr& testTP =
0111 ttClusterAssociationMapHandle->findTrackingParticlePtr(tempStubRef->clusterRef(0));
0112
0113
0114
0115 stubToTrackingParticleMap.emplace(tempStubRef, testTP);
0116
0117
0118 continue;
0119 } else {
0120
0121 continue;
0122 }
0123 }
0124 else {
0125
0126 TrackingParticle* prevTPAddress = nullptr;
0127 unsigned int whichTP = 0;
0128
0129 const std::vector<TrackingParticlePtr>& trackingParticles0 =
0130 ttClusterAssociationMapHandle->findTrackingParticlePtrs(tempStubRef->clusterRef(0));
0131 const std::vector<TrackingParticlePtr>& trackingParticles1 =
0132 ttClusterAssociationMapHandle->findTrackingParticlePtrs(tempStubRef->clusterRef(1));
0133
0134 bool escape = false;
0135
0136 for (unsigned int i = 0; i < trackingParticles0.size() && !escape; i++) {
0137
0138 if (trackingParticles0.at(i).isNull())
0139 continue;
0140
0141 for (unsigned int k = 0; k < trackingParticles1.size() && !escape; k++) {
0142
0143 if (trackingParticles1.at(k).isNull())
0144 continue;
0145
0146 if (trackingParticles0.at(i).get() == trackingParticles1.at(k).get()) {
0147
0148 if (prevTPAddress == nullptr) {
0149 prevTPAddress = const_cast<TrackingParticle*>(trackingParticles1.at(k).get());
0150 whichTP = k;
0151 }
0152
0153
0154
0155 if (prevTPAddress != const_cast<TrackingParticle*>(trackingParticles1.at(k).get())) {
0156 escape = true;
0157 continue;
0158 }
0159 }
0160 }
0161 }
0162
0163
0164
0165 if (escape)
0166 continue;
0167
0168 if (prevTPAddress == nullptr) {
0169
0170 continue;
0171 } else {
0172
0173
0174
0175
0176 TrackingParticlePtr testTP = trackingParticles1.at(whichTP);
0177
0178
0179
0180 stubToTrackingParticleMap.emplace(tempStubRef, testTP);
0181
0182
0183 continue;
0184 }
0185 }
0186 }
0187 }
0188 }
0189 }
0190
0191
0192
0193 for (auto& p : trackingParticleToStubVectorMap) {
0194
0195 std::vector<TTStubRef>& tempVector = p.second;
0196
0197
0198 std::sort(tempVector.begin(), tempVector.end());
0199 tempVector.erase(std::unique(tempVector.begin(), tempVector.end()), tempVector.end());
0200 }
0201
0202
0203 edm::RefProd<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> theCluAssoMap(ttClusterAssociationMapHandle);
0204
0205
0206 associationMapForOutput->setTTStubToTrackingParticleMap(stubToTrackingParticleMap);
0207 associationMapForOutput->setTrackingParticleToTTStubsMap(trackingParticleToStubVectorMap);
0208 associationMapForOutput->setTTClusterAssociationMap(theCluAssoMap);
0209
0210
0211 iEvent.put(std::move(associationMapForOutput), ttStubsInputTags_.at(ncont1).instance());
0212
0213 ++ncont1;
0214 }
0215 }