Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-12 23:42:06

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