Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:16

0001 #include "FWCore/Framework/interface/global/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/ESHandle.h"
0007 #include "FWCore/Framework/interface/ConsumesCollector.h"
0008 #include "FWCore/PluginManager/interface/ModuleDef.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/Framework/interface/ESHandle.h"
0012 
0013 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0014 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0015 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0016 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0017 
0018 #include "DataFormats/Common/interface/DetSetVector.h"
0019 #include "DataFormats/Phase2TrackerCluster/interface/Phase2TrackerCluster1D.h"
0020 #include "DataFormats/TrackerRecHit2D/interface/Phase2TrackerRecHit1D.h"
0021 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0022 #include "DataFormats/DetId/interface/DetId.h"
0023 
0024 #include "RecoLocalTracker/Records/interface/TkPhase2OTCPERecord.h"
0025 #include "RecoLocalTracker/Phase2TrackerRecHits/interface/Phase2StripCPE.h"
0026 
0027 #include <vector>
0028 #include <string>
0029 
0030 class Phase2TrackerRecHits : public edm::global::EDProducer<> {
0031 public:
0032   explicit Phase2TrackerRecHits(const edm::ParameterSet& conf);
0033   ~Phase2TrackerRecHits() override{};
0034   void produce(edm::StreamID sid, edm::Event& event, const edm::EventSetup& eventSetup) const final;
0035 
0036 private:
0037   edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> const tTrackerGeom_;
0038   edm::ESGetToken<ClusterParameterEstimator<Phase2TrackerCluster1D>, TkPhase2OTCPERecord> const tCPE_;
0039 
0040   edm::EDGetTokenT<Phase2TrackerCluster1DCollectionNew> token_;
0041 };
0042 
0043 Phase2TrackerRecHits::Phase2TrackerRecHits(edm::ParameterSet const& conf)
0044     : tTrackerGeom_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0045       tCPE_(esConsumes(conf.getParameter<edm::ESInputTag>("Phase2StripCPE"))),
0046       token_(consumes<Phase2TrackerCluster1DCollectionNew>(conf.getParameter<edm::InputTag>("src"))) {
0047   produces<Phase2TrackerRecHit1DCollectionNew>();
0048 }
0049 
0050 void Phase2TrackerRecHits::produce(edm::StreamID sid, edm::Event& event, const edm::EventSetup& eventSetup) const {
0051   // Get the Clusters
0052   edm::Handle<Phase2TrackerCluster1DCollectionNew> clusters;
0053   event.getByToken(token_, clusters);
0054 
0055   // load the cpe via the eventsetup
0056   const auto& cpe = &eventSetup.getData(tCPE_);
0057 
0058   // Get the geometry
0059   const TrackerGeometry* tkGeom = &eventSetup.getData(tTrackerGeom_);
0060 
0061   // Global container for the RecHits of each module
0062   auto outputRecHits = std::make_unique<Phase2TrackerRecHit1DCollectionNew>();
0063 
0064   // Loop over clusters
0065   for (const auto& clusterDetSet : *clusters) {
0066     DetId detId(clusterDetSet.detId());
0067 
0068     // Geometry
0069     const GeomDetUnit* geomDetUnit(tkGeom->idToDetUnit(detId));
0070 
0071     // Container for the clusters that will be produced for this modules
0072     Phase2TrackerRecHit1DCollectionNew::FastFiller rechits(*outputRecHits, clusterDetSet.detId());
0073 
0074     for (const auto& clusterRef : clusterDetSet) {
0075       ClusterParameterEstimator<Phase2TrackerCluster1D>::LocalValues lv =
0076           cpe->localParameters(clusterRef, *geomDetUnit);
0077 
0078       // Create a persistent edm::Ref to the cluster
0079       edm::Ref<Phase2TrackerCluster1DCollectionNew, Phase2TrackerCluster1D> cluster =
0080           edmNew::makeRefTo(clusters, &clusterRef);
0081 
0082       // Make a RecHit and add it to the DetSet
0083       Phase2TrackerRecHit1D hit(lv.first, lv.second, *geomDetUnit, cluster);
0084 
0085       rechits.push_back(hit);
0086     }
0087   }
0088 
0089   outputRecHits->shrink_to_fit();
0090   event.put(std::move(outputRecHits));
0091 }
0092 
0093 DEFINE_FWK_MODULE(Phase2TrackerRecHits);