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/SiStripCluster/interface/SiStripClusterTools.h"
0008 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
0009 
0010 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0011 
0012 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
0013 
0014 #include "RecoTracker/MkFit/interface/MkFitHitWrapper.h"
0015 #include "RecoTracker/MkFit/interface/MkFitClusterIndexToHit.h"
0016 #include "RecoTracker/MkFit/interface/MkFitGeometry.h"
0017 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
0018 
0019 #include "convertHits.h"
0020 
0021 namespace {
0022   class ConvertHitTraits {
0023   public:
0024     ConvertHitTraits(float minCharge) : minGoodStripCharge_(minCharge) {}
0025     using Clusters = edmNew::DetSetVector<SiStripCluster>;
0026     using Cluster = Clusters::data_type;
0027 
0028     static constexpr bool applyCCC() { return true; }
0029     static float chargeScale(DetId id) { return siStripClusterTools::sensorThicknessInverse(id); }
0030     static const Cluster& cluster(const Clusters& prod, unsigned int index) { return prod.data()[index]; }
0031     static float clusterCharge(const Cluster& clu, float scale) { return clu.charge() * scale; }
0032     bool passCCC(float charge) const { return charge > minGoodStripCharge_; }
0033     static void setDetails(mkfit::Hit& mhit, const Cluster& clu, int shortId, float charge) {
0034       mhit.setupAsStrip(shortId, charge, clu.size());
0035     }
0036 
0037   private:
0038     const float minGoodStripCharge_;
0039   };
0040 }  // namespace
0041 
0042 class MkFitSiStripHitConverter : public edm::global::EDProducer<> {
0043 public:
0044   explicit MkFitSiStripHitConverter(edm::ParameterSet const& iConfig);
0045   ~MkFitSiStripHitConverter() override = default;
0046 
0047   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0048 
0049 private:
0050   void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0051 
0052   const edm::EDGetTokenT<SiStripRecHit2DCollection> stripRphiRecHitToken_;
0053   const edm::EDGetTokenT<SiStripRecHit2DCollection> stripStereoRecHitToken_;
0054   const edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster>> stripClusterToken_;
0055   const edm::ESGetToken<TransientTrackingRecHitBuilder, TransientRecHitRecord> ttrhBuilderToken_;
0056   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> ttopoToken_;
0057   const edm::ESGetToken<MkFitGeometry, TrackerRecoGeometryRecord> mkFitGeomToken_;
0058   const edm::EDPutTokenT<MkFitHitWrapper> wrapperPutToken_;
0059   const edm::EDPutTokenT<MkFitClusterIndexToHit> clusterIndexPutToken_;
0060   const edm::EDPutTokenT<std::vector<int>> layerIndexPutToken_;
0061   const edm::EDPutTokenT<std::vector<float>> clusterChargePutToken_;
0062   const ConvertHitTraits convertTraits_;
0063 };
0064 
0065 MkFitSiStripHitConverter::MkFitSiStripHitConverter(edm::ParameterSet const& iConfig)
0066     : stripRphiRecHitToken_{consumes(iConfig.getParameter<edm::InputTag>("rphiHits"))},
0067       stripStereoRecHitToken_{consumes(iConfig.getParameter<edm::InputTag>("stereoHits"))},
0068       stripClusterToken_{consumes(iConfig.getParameter<edm::InputTag>("clusters"))},
0069       ttrhBuilderToken_{esConsumes(iConfig.getParameter<edm::ESInputTag>("ttrhBuilder"))},
0070       ttopoToken_{esConsumes()},
0071       mkFitGeomToken_{esConsumes()},
0072       wrapperPutToken_{produces()},
0073       clusterIndexPutToken_{produces()},
0074       layerIndexPutToken_{produces()},
0075       clusterChargePutToken_{produces()},
0076       convertTraits_{static_cast<float>(
0077           iConfig.getParameter<edm::ParameterSet>("minGoodStripCharge").getParameter<double>("value"))} {}
0078 
0079 void MkFitSiStripHitConverter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0080   edm::ParameterSetDescription desc;
0081 
0082   desc.add("rphiHits", edm::InputTag{"siStripMatchedRecHits", "rphiRecHit"});
0083   desc.add("stereoHits", edm::InputTag{"siStripMatchedRecHits", "stereoRecHit"});
0084   desc.add("clusters", edm::InputTag{"siStripClusters"});
0085   desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
0086 
0087   edm::ParameterSetDescription descCCC;
0088   descCCC.add<double>("value");
0089   desc.add("minGoodStripCharge", descCCC);
0090 
0091   descriptions.add("mkFitSiStripHitConverterDefault", desc);
0092 }
0093 
0094 void MkFitSiStripHitConverter::produce(edm::StreamID iID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0095   const auto& ttrhBuilder = iSetup.getData(ttrhBuilderToken_);
0096   const auto& ttopo = iSetup.getData(ttopoToken_);
0097   const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
0098 
0099   MkFitHitWrapper hitWrapper;
0100   MkFitClusterIndexToHit clusterIndexToHit;
0101   std::vector<int> layerIndexToHit;
0102   std::vector<float> clusterCharge;
0103 
0104   edm::ProductID stripClusterID;
0105   const auto& stripRphiHits = iEvent.get(stripRphiRecHitToken_);
0106   const auto& stripStereoHits = iEvent.get(stripStereoRecHitToken_);
0107   const auto maxSizeGuess(stripRphiHits.dataSize() + stripStereoHits.dataSize());
0108   auto const& clusters = iEvent.get(stripClusterToken_);
0109 
0110   auto convert = [&](auto& hits) {
0111     return mkfit::convertHits(convertTraits_,
0112                               hits,
0113                               clusters,
0114                               hitWrapper.hits(),
0115                               clusterIndexToHit.hits(),
0116                               layerIndexToHit,
0117                               clusterCharge,
0118                               ttopo,
0119                               ttrhBuilder,
0120                               mkFitGeom,
0121                               maxSizeGuess);
0122   };
0123 
0124   if (not stripRphiHits.empty()) {
0125     stripClusterID = convert(stripRphiHits);
0126   }
0127   if (not stripStereoHits.empty()) {
0128     auto stripStereoClusterID = convert(stripStereoHits);
0129     if (stripRphiHits.empty()) {
0130       stripClusterID = stripStereoClusterID;
0131     } else if (stripClusterID != stripStereoClusterID) {
0132       throw cms::Exception("LogicError") << "Encountered different cluster ProductIDs for strip RPhi hits ("
0133                                          << stripClusterID << ") and stereo (" << stripStereoClusterID << ")";
0134     }
0135   }
0136 
0137   hitWrapper.setClustersID(stripClusterID);
0138 
0139   iEvent.emplace(wrapperPutToken_, std::move(hitWrapper));
0140   iEvent.emplace(clusterIndexPutToken_, std::move(clusterIndexToHit));
0141   iEvent.emplace(layerIndexPutToken_, std::move(layerIndexToHit));
0142   iEvent.emplace(clusterChargePutToken_, std::move(clusterCharge));
0143 }
0144 
0145 DEFINE_FWK_MODULE(MkFitSiStripHitConverter);