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