File indexing completed on 2023-03-17 11:19:31
0001 #include "FWCore/Framework/interface/stream/EDProducer.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/EventSetup.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/Utilities/interface/InputTag.h"
0006 #include "FWCore/Framework/interface/ConsumesCollector.h"
0007 #include "FWCore/PluginManager/interface/ModuleDef.h"
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011
0012 #ifdef VERIFY_PH2_TK_CLUS
0013 #include "Phase2TrackerClusterizerAlgorithm.h"
0014 #endif
0015 #include "Phase2TrackerClusterizerSequentialAlgorithm.h"
0016
0017 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0018 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0019 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0020 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0021
0022 #include "DataFormats/Common/interface/DetSetVector.h"
0023 #include "DataFormats/Phase2TrackerDigi/interface/Phase2TrackerDigi.h"
0024 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0025 #include "DataFormats/DetId/interface/DetId.h"
0026
0027 #include <vector>
0028 #include <memory>
0029
0030 class Phase2TrackerClusterizer : public edm::stream::EDProducer<> {
0031 public:
0032 explicit Phase2TrackerClusterizer(const edm::ParameterSet& conf);
0033 ~Phase2TrackerClusterizer() override;
0034 void produce(edm::Event& event, const edm::EventSetup& eventSetup) override;
0035
0036 private:
0037 #ifdef VERIFY_PH2_TK_CLUS
0038 std::unique_ptr<Phase2TrackerClusterizerAlgorithm> clusterizer_;
0039 #endif
0040 edm::EDGetTokenT<edm::DetSetVector<Phase2TrackerDigi> > token_;
0041 };
0042
0043
0044
0045
0046
0047 Phase2TrackerClusterizer::Phase2TrackerClusterizer(edm::ParameterSet const& conf)
0048 :
0049 #ifdef VERIFY_PH2_TK_CLUS
0050 clusterizer_(new Phase2TrackerClusterizerAlgorithm(conf.getParameter<unsigned int>("maxClusterSize"),
0051 conf.getParameter<unsigned int>("maxNumberClusters"))),
0052 #endif
0053 token_(consumes<edm::DetSetVector<Phase2TrackerDigi> >(conf.getParameter<edm::InputTag>("src"))) {
0054 produces<Phase2TrackerCluster1DCollectionNew>();
0055 }
0056
0057 Phase2TrackerClusterizer::~Phase2TrackerClusterizer() {}
0058
0059
0060
0061
0062
0063 void Phase2TrackerClusterizer::produce(edm::Event& event, const edm::EventSetup& eventSetup) {
0064
0065 edm::Handle<edm::DetSetVector<Phase2TrackerDigi> > digis;
0066 event.getByToken(token_, digis);
0067
0068 #ifdef VERIFY_PH2_TK_CLUS
0069
0070 edm::ESHandle<TrackerGeometry> geomHandle;
0071 eventSetup.get<TrackerDigiGeometryRecord>().get(geomHandle);
0072 const TrackerGeometry* tkGeom(&(*geomHandle));
0073
0074 auto outputClustersOld = std::make_unique<Phase2TrackerCluster1DCollectionNew>();
0075 #endif
0076 auto outputClusters = std::make_unique<Phase2TrackerCluster1DCollectionNew>();
0077
0078
0079 for (const auto& DSViter : *digis) {
0080 DetId detId(DSViter.detId());
0081
0082 Phase2TrackerCluster1DCollectionNew::FastFiller clusters(*outputClusters, DSViter.detId());
0083 Phase2TrackerClusterizerSequentialAlgorithm algo;
0084 algo.clusterizeDetUnit(DSViter, clusters);
0085 if (clusters.empty())
0086 clusters.abort();
0087
0088 #ifdef VERIFY_PH2_TK_CLUS
0089 if (!clusters.empty()) {
0090 auto cp = clusters[0].column();
0091 auto sp = clusters[0].firstStrip();
0092 for (auto const& cl : clusters) {
0093 if (cl.column() < cp)
0094 std::cout << "column not in order! " << std::endl;
0095 if (cl.column() == cp && cl.firstStrip() < sp)
0096 std::cout << "strip not in order! " << std::endl;
0097 cp = cl.column();
0098 sp = cl.firstStrip();
0099 }
0100 }
0101 #endif
0102
0103 #ifdef VERIFY_PH2_TK_CLUS
0104
0105 const GeomDetUnit* geomDetUnit(tkGeom->idToDetUnit(detId));
0106 const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(geomDetUnit);
0107 if (!pixDet)
0108 assert(0);
0109
0110
0111 Phase2TrackerCluster1DCollectionNew::FastFiller clustersOld(*outputClustersOld, DSViter.detId());
0112
0113
0114 clusterizer_->setup(pixDet);
0115
0116
0117
0118 clusterizer_->clusterizeDetUnit(DSViter, clustersOld);
0119 if (clustersOld.empty())
0120 clustersOld.abort();
0121
0122 if (clusters.size() != clustersOld.size()) {
0123 std::cout << "SIZEs " << int(detId) << ' ' << clusters.size() << ' ' << clustersOld.size() << std::endl;
0124 for (auto const& cl : clusters)
0125 std::cout << cl.size() << ' ' << cl.threshold() << ' ' << cl.firstRow() << ' ' << cl.column() << std::endl;
0126 std::cout << "Old " << std::endl;
0127 for (auto const& cl : clustersOld)
0128 std::cout << cl.size() << ' ' << cl.threshold() << ' ' << cl.firstRow() << ' ' << cl.column() << std::endl;
0129 }
0130 #endif
0131 }
0132
0133 #ifdef VERIFY_PH2_TK_CLUS
0134
0135 assert(outputClusters->dataSize() == outputClustersOld->dataSize());
0136 for (auto i = 0U; i < outputClusters->dataSize(); ++i) {
0137 assert(outputClusters->data()[i].size() == outputClustersOld->data()[i].size());
0138 assert(outputClusters->data()[i].threshold() == outputClustersOld->data()[i].threshold());
0139 assert(outputClusters->data()[i].firstRow() == outputClustersOld->data()[i].firstRow());
0140 assert(outputClusters->data()[i].column() == outputClustersOld->data()[i].column());
0141 }
0142 #endif
0143
0144
0145 outputClusters->shrink_to_fit();
0146 event.put(std::move(outputClusters));
0147 }
0148
0149 DEFINE_FWK_MODULE(Phase2TrackerClusterizer);