File indexing completed on 2024-04-06 12:01:10
0001 #include "CommonTools/RecoAlgos/interface/ClusterStorer.h"
0002
0003 #include "FWCore/Utilities/interface/Exception.h"
0004
0005 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0006 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0007 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
0008 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0009 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
0010 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0011 #include "DataFormats/TrackerRecHit2D/interface/Phase2TrackerRecHit1D.h"
0012 #include "DataFormats/TrackerRecHit2D/interface/VectorHit.h"
0013
0014 #include "DataFormats/TrackerRecHit2D/interface/FastTrackerRecHit.h"
0015 #include "DataFormats/TrackerRecHit2D/interface/FastProjectedTrackerRecHit.h"
0016 #include "DataFormats/TrackerRecHit2D/interface/FastMatchedTrackerRecHit.h"
0017
0018 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
0019
0020 namespace helper {
0021
0022
0023
0024 void ClusterStorer::addCluster(TrackingRecHitCollection &hits, size_t index) {
0025 TrackingRecHit &newHit = hits[index];
0026 const std::type_info &hit_type = typeid(newHit);
0027 if (hit_type == typeid(SiPixelRecHit)) {
0028
0029 pixelClusterRecords_.push_back(PixelClusterHitRecord(static_cast<SiPixelRecHit &>(newHit), hits, index));
0030 } else if (hit_type == typeid(SiStripRecHit1D)) {
0031
0032 stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit1D &>(newHit), hits, index));
0033 } else if (hit_type == typeid(SiStripRecHit2D)) {
0034
0035 stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit2D &>(newHit), hits, index));
0036 } else if (hit_type == typeid(SiStripMatchedRecHit2D)) {
0037
0038 SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D &>(newHit);
0039 stripClusterRecords_.push_back(StripClusterHitRecord(mhit.monoHit(), hits, index));
0040 stripClusterRecords_.push_back(StripClusterHitRecord(mhit.stereoHit(), hits, index));
0041 } else if (hit_type == typeid(ProjectedSiStripRecHit2D)) {
0042
0043 ProjectedSiStripRecHit2D &phit = static_cast<ProjectedSiStripRecHit2D &>(newHit);
0044 stripClusterRecords_.push_back(StripClusterHitRecord(phit.originalHit(), hits, index));
0045 } else if (hit_type == typeid(Phase2TrackerRecHit1D)) {
0046
0047
0048 phase2OTClusterRecords_.push_back(
0049 Phase2OTClusterHitRecord(static_cast<Phase2TrackerRecHit1D &>(newHit), hits, index));
0050 } else if (hit_type == typeid(VectorHit)) {
0051
0052
0053 phase2OTClusterRecords_.push_back(Phase2OTClusterHitRecord(static_cast<VectorHit &>(newHit), hits, index));
0054 } else {
0055 if (hit_type == typeid(FastTrackerRecHit) || hit_type == typeid(FastProjectedTrackerRecHit) ||
0056 hit_type == typeid(FastMatchedTrackerRecHit)) {
0057
0058
0059
0060
0061
0062 } else {
0063
0064 throw cms::Exception("UnknownHitType") << "helper::ClusterStorer::addCluster: "
0065 << "Unknown hit type " << hit_type.name() << ".\n";
0066 }
0067 }
0068 }
0069
0070
0071 void ClusterStorer::clear() {
0072 pixelClusterRecords_.clear();
0073 stripClusterRecords_.clear();
0074 phase2OTClusterRecords_.clear();
0075 }
0076
0077
0078 void ClusterStorer::processAllClusters(
0079 edmNew::DetSetVector<SiPixelCluster> &pixelDsvToFill,
0080 edm::RefProd<edmNew::DetSetVector<SiPixelCluster> > refPixelClusters,
0081 edmNew::DetSetVector<SiStripCluster> &stripDsvToFill,
0082 edm::RefProd<edmNew::DetSetVector<SiStripCluster> > refStripClusters,
0083 edmNew::DetSetVector<Phase2TrackerCluster1D> &phase2OTDsvToFill,
0084 edm::RefProd<edmNew::DetSetVector<Phase2TrackerCluster1D> > refPhase2OTClusters) {
0085 if (!pixelClusterRecords_.empty()) {
0086 this->processClusters<SiPixelRecHit, SiPixelCluster>(pixelClusterRecords_, pixelDsvToFill, refPixelClusters);
0087 }
0088 if (!stripClusterRecords_.empty()) {
0089
0090
0091
0092
0093
0094 this->processClusters<SiStripRecHit2D, SiStripCluster>(stripClusterRecords_, stripDsvToFill, refStripClusters);
0095 }
0096 if (!phase2OTClusterRecords_.empty()) {
0097 this->processClusters<Phase2TrackerRecHit1D, Phase2TrackerCluster1D>(
0098 phase2OTClusterRecords_, phase2OTDsvToFill, refPhase2OTClusters);
0099 }
0100 }
0101
0102
0103 template <typename HitType, typename ClusterType>
0104 void ClusterStorer::processClusters(std::vector<ClusterHitRecord<typename HitType::ClusterRef> > &clusterRecords,
0105 edmNew::DetSetVector<ClusterType> &dsvToFill,
0106 edm::RefProd<edmNew::DetSetVector<ClusterType> > &refprod) {
0107 std::sort(clusterRecords.begin(), clusterRecords.end());
0108 typedef typename std::vector<ClusterHitRecord<typename HitType::ClusterRef> >::iterator RIT;
0109 RIT it = clusterRecords.begin(), end = clusterRecords.end();
0110 size_t clusters = 0;
0111 while (it != end) {
0112 RIT it2 = it;
0113 uint32_t detid = it->detid();
0114
0115
0116 while ((it2 != end) && (it2->detid() == detid)) {
0117 ++it2;
0118 }
0119
0120
0121
0122 typename edmNew::DetSetVector<ClusterType>::FastFiller filler(dsvToFill, detid);
0123 typename HitType::ClusterRef lastRef, newRef;
0124 for (; it != it2; ++it) {
0125
0126 if (it->clusterRef() != lastRef) {
0127 lastRef = it->clusterRef();
0128
0129 filler.push_back(*lastRef);
0130
0131 newRef = typename HitType::ClusterRef(refprod, clusters++);
0132 }
0133 it->template rekey<HitType>(newRef);
0134 }
0135
0136 }
0137 }
0138
0139
0140
0141
0142
0143
0144 template <typename ClusterRefType>
0145 template <typename RecHitType>
0146 void ClusterStorer::ClusterHitRecord<ClusterRefType>::rekey(const ClusterRefType &newRef) {
0147 TrackingRecHit &genericHit = (*hits_)[index_];
0148 RecHitType *hit = nullptr;
0149 if (genericHit.geographicalId().rawId() == detid_) {
0150 hit = dynamic_cast<RecHitType *>(&genericHit);
0151 }
0152 assert(hit != nullptr);
0153 assert(hit->cluster() == ref_);
0154 hit->setClusterRef(newRef);
0155 }
0156
0157
0158
0159
0160 template <>
0161 template <typename RecHitType>
0162 void ClusterStorer::ClusterHitRecord<SiStripRecHit2D::ClusterRef>::rekey(const SiStripRecHit2D::ClusterRef &newRef) {
0163 TrackingRecHit &genericHit = (*hits_)[index_];
0164 const std::type_info &hit_type = typeid(genericHit);
0165
0166 OmniClusterRef *cluRef = nullptr;
0167 if (typeid(SiStripRecHit1D) == hit_type) {
0168 cluRef = &static_cast<SiStripRecHit1D &>(genericHit).omniCluster();
0169 } else if (typeid(SiStripRecHit2D) == hit_type) {
0170 cluRef = &static_cast<SiStripRecHit2D &>(genericHit).omniCluster();
0171 } else if (typeid(SiStripMatchedRecHit2D) == hit_type) {
0172 SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D &>(genericHit);
0173 cluRef = (SiStripDetId(detid_).stereo() ? &mhit.stereoClusterRef() : &mhit.monoClusterRef());
0174 } else if (typeid(ProjectedSiStripRecHit2D) == hit_type) {
0175 cluRef = &static_cast<ProjectedSiStripRecHit2D &>(genericHit).omniCluster();
0176 }
0177
0178 assert(cluRef != nullptr);
0179 assert(cluRef->key() == ref_.key());
0180 (*cluRef) = OmniClusterRef(newRef);
0181 }
0182
0183
0184
0185 template <>
0186 template <typename RecHitType>
0187 void ClusterStorer::ClusterHitRecord<Phase2TrackerRecHit1D::ClusterRef>::rekey(
0188 const Phase2TrackerRecHit1D::ClusterRef &newRef) {
0189 TrackingRecHit &genericHit = (*hits_)[index_];
0190 const std::type_info &hit_type = typeid(genericHit);
0191
0192 OmniClusterRef *cluRef = nullptr;
0193 if (typeid(Phase2TrackerRecHit1D) == hit_type) {
0194 cluRef = &static_cast<Phase2TrackerRecHit1D &>(genericHit).omniCluster();
0195 } else if (typeid(VectorHit) == hit_type) {
0196 VectorHit &vHit = static_cast<VectorHit &>(genericHit);
0197
0198
0199 cluRef = (SiStripDetId(detid_).stereo() ? &vHit.upperClusterRef() : &vHit.lowerClusterRef());
0200 }
0201
0202 assert(cluRef != nullptr);
0203 assert(cluRef->key() == ref_.key());
0204 (*cluRef) = OmniClusterRef(newRef);
0205 }
0206
0207 }