Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:13

0001 #include "FWCore/Framework/interface/global/EDProducer.h"
0002 
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 
0007 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0008 
0009 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0010 
0011 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
0012 
0013 #include "RecoTracker/MkFit/interface/MkFitHitWrapper.h"
0014 #include "RecoTracker/MkFit/interface/MkFitClusterIndexToHit.h"
0015 #include "RecoTracker/MkFit/interface/MkFitGeometry.h"
0016 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
0017 
0018 #include "convertHits.h"
0019 
0020 namespace {
0021   struct ConvertHitTraits {
0022     static constexpr bool applyCCC() { return false; }
0023     static std::nullptr_t clusterCharge(const SiPixelRecHit& hit, DetId hitId) { return nullptr; }
0024     static bool passCCC(std::nullptr_t) { return true; }
0025     static void setDetails(mkfit::Hit& mhit, const SiPixelCluster& cluster, int shortId, std::nullptr_t) {
0026       mhit.setupAsPixel(shortId, cluster.sizeX(), cluster.sizeY());
0027     }
0028   };
0029 }  // namespace
0030 
0031 class MkFitSiPixelHitConverter : public edm::global::EDProducer<> {
0032 public:
0033   explicit MkFitSiPixelHitConverter(edm::ParameterSet const& iConfig);
0034   ~MkFitSiPixelHitConverter() override = default;
0035 
0036   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0037 
0038 private:
0039   void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0040 
0041   const edm::EDGetTokenT<SiPixelRecHitCollection> pixelRecHitToken_;
0042   const edm::ESGetToken<TransientTrackingRecHitBuilder, TransientRecHitRecord> ttrhBuilderToken_;
0043   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> ttopoToken_;
0044   const edm::ESGetToken<MkFitGeometry, TrackerRecoGeometryRecord> mkFitGeomToken_;
0045   const edm::EDPutTokenT<MkFitHitWrapper> wrapperPutToken_;
0046   const edm::EDPutTokenT<MkFitClusterIndexToHit> clusterIndexPutToken_;
0047 };
0048 
0049 MkFitSiPixelHitConverter::MkFitSiPixelHitConverter(edm::ParameterSet const& iConfig)
0050     : pixelRecHitToken_{consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("hits"))},
0051       ttrhBuilderToken_{esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(
0052           iConfig.getParameter<edm::ESInputTag>("ttrhBuilder"))},
0053       ttopoToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
0054       mkFitGeomToken_{esConsumes<MkFitGeometry, TrackerRecoGeometryRecord>()},
0055       wrapperPutToken_{produces<MkFitHitWrapper>()},
0056       clusterIndexPutToken_{produces<MkFitClusterIndexToHit>()} {}
0057 
0058 void MkFitSiPixelHitConverter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0059   edm::ParameterSetDescription desc;
0060 
0061   desc.add("hits", edm::InputTag{"siPixelRecHits"});
0062   desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
0063 
0064   descriptions.addWithDefaultLabel(desc);
0065 }
0066 
0067 void MkFitSiPixelHitConverter::produce(edm::StreamID iID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0068   const auto& ttrhBuilder = iSetup.getData(ttrhBuilderToken_);
0069   const auto& ttopo = iSetup.getData(ttopoToken_);
0070   const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
0071 
0072   MkFitHitWrapper hitWrapper;
0073   MkFitClusterIndexToHit clusterIndexToHit;
0074 
0075   std::vector<float> dummy;
0076   auto pixelClusterID = mkfit::convertHits(ConvertHitTraits{},
0077                                            iEvent.get(pixelRecHitToken_),
0078                                            hitWrapper.hits(),
0079                                            clusterIndexToHit.hits(),
0080                                            dummy,
0081                                            ttopo,
0082                                            ttrhBuilder,
0083                                            mkFitGeom);
0084 
0085   hitWrapper.setClustersID(pixelClusterID);
0086 
0087   iEvent.emplace(wrapperPutToken_, std::move(hitWrapper));
0088   iEvent.emplace(clusterIndexPutToken_, std::move(clusterIndexToHit));
0089 }
0090 
0091 DEFINE_FWK_MODULE(MkFitSiPixelHitConverter);