Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-18 03:42:17

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     using Clusters = SiPixelClusterCollectionNew;
0023     using Cluster = Clusters::data_type;
0024 
0025     static constexpr bool applyCCC() { return false; }
0026     static float chargeScale(DetId) { return 0; }
0027     static const Cluster& cluster(const Clusters& prod, unsigned int index) { return prod.data()[index]; }
0028     static std::nullptr_t clusterCharge(const Cluster&, float) { return nullptr; }
0029     static bool passCCC(std::nullptr_t) { return true; }
0030     static void setDetails(mkfit::Hit& mhit, const Cluster& clu, int shortId, std::nullptr_t) {
0031       mhit.setupAsPixel(shortId, clu.sizeX(), clu.sizeY());
0032     }
0033   };
0034 }  // namespace
0035 
0036 class MkFitSiPixelHitConverter : public edm::global::EDProducer<> {
0037 public:
0038   explicit MkFitSiPixelHitConverter(edm::ParameterSet const& iConfig);
0039   ~MkFitSiPixelHitConverter() override = default;
0040 
0041   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0042 
0043 private:
0044   void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0045 
0046   const edm::EDGetTokenT<SiPixelRecHitCollection> pixelRecHitToken_;
0047   const edm::EDGetTokenT<SiPixelClusterCollectionNew> pixelClusterToken_;
0048   const edm::ESGetToken<TransientTrackingRecHitBuilder, TransientRecHitRecord> ttrhBuilderToken_;
0049   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> ttopoToken_;
0050   const edm::ESGetToken<MkFitGeometry, TrackerRecoGeometryRecord> mkFitGeomToken_;
0051   const edm::EDPutTokenT<MkFitHitWrapper> wrapperPutToken_;
0052   const edm::EDPutTokenT<MkFitClusterIndexToHit> clusterIndexPutToken_;
0053   const edm::EDPutTokenT<std::vector<int>> layerIndexPutToken_;
0054 };
0055 
0056 MkFitSiPixelHitConverter::MkFitSiPixelHitConverter(edm::ParameterSet const& iConfig)
0057     : pixelRecHitToken_{consumes(iConfig.getParameter<edm::InputTag>("hits"))},
0058       pixelClusterToken_{consumes(iConfig.getParameter<edm::InputTag>("clusters"))},
0059       ttrhBuilderToken_{esConsumes(iConfig.getParameter<edm::ESInputTag>("ttrhBuilder"))},
0060       ttopoToken_{esConsumes()},
0061       mkFitGeomToken_{esConsumes()},
0062       wrapperPutToken_{produces()},
0063       clusterIndexPutToken_{produces()},
0064       layerIndexPutToken_{produces()} {}
0065 
0066 void MkFitSiPixelHitConverter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0067   edm::ParameterSetDescription desc;
0068 
0069   desc.add("hits", edm::InputTag{"siPixelRecHits"});
0070   desc.add("clusters", edm::InputTag{"siPixelClusters"});
0071   desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
0072 
0073   descriptions.addWithDefaultLabel(desc);
0074 }
0075 
0076 void MkFitSiPixelHitConverter::produce(edm::StreamID iID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0077   const auto& ttrhBuilder = iSetup.getData(ttrhBuilderToken_);
0078   const auto& ttopo = iSetup.getData(ttopoToken_);
0079   const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
0080 
0081   MkFitHitWrapper hitWrapper;
0082   MkFitClusterIndexToHit clusterIndexToHit;
0083   std::vector<int> layerIndexToHit;
0084 
0085   std::vector<float> dummy;
0086   auto const& hits = iEvent.get(pixelRecHitToken_);
0087   auto pixelClusterID = mkfit::convertHits(ConvertHitTraits{},
0088                                            hits,
0089                                            iEvent.get(pixelClusterToken_),
0090                                            hitWrapper.hits(),
0091                                            clusterIndexToHit.hits(),
0092                                            layerIndexToHit,
0093                                            dummy,
0094                                            ttopo,
0095                                            ttrhBuilder,
0096                                            mkFitGeom,
0097                                            hits.dataSize());
0098 
0099   hitWrapper.setClustersID(pixelClusterID);
0100 
0101   iEvent.emplace(wrapperPutToken_, std::move(hitWrapper));
0102   iEvent.emplace(clusterIndexPutToken_, std::move(clusterIndexToHit));
0103   iEvent.emplace(layerIndexPutToken_, std::move(layerIndexToHit));
0104 }
0105 
0106 DEFINE_FWK_MODULE(MkFitSiPixelHitConverter);