Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-25 02:14:11

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