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/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
0026 static constexpr bool applyCCC() { return true; }
0027 static float clusterCharge(const SiStripRecHit2D& hit, DetId hitId) {
0028 return siStripClusterTools::chargePerCM(hitId, hit.firstClusterRef().stripCluster());
0029 }
0030 bool passCCC(float charge) const { return charge > minGoodStripCharge_; }
0031 static void setDetails(mkfit::Hit& mhit, const SiStripCluster& cluster, int shortId, float charge) {
0032 mhit.setupAsStrip(shortId, charge, cluster.amplitudes().size());
0033 }
0034
0035 private:
0036 const float minGoodStripCharge_;
0037 };
0038 }
0039
0040 class MkFitSiStripHitConverter : public edm::global::EDProducer<> {
0041 public:
0042 explicit MkFitSiStripHitConverter(edm::ParameterSet const& iConfig);
0043 ~MkFitSiStripHitConverter() override = default;
0044
0045 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0046
0047 private:
0048 void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0049
0050 const edm::EDGetTokenT<SiStripRecHit2DCollection> stripRphiRecHitToken_;
0051 const edm::EDGetTokenT<SiStripRecHit2DCollection> stripStereoRecHitToken_;
0052 const edm::ESGetToken<TransientTrackingRecHitBuilder, TransientRecHitRecord> ttrhBuilderToken_;
0053 const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> ttopoToken_;
0054 const edm::ESGetToken<MkFitGeometry, TrackerRecoGeometryRecord> mkFitGeomToken_;
0055 const edm::EDPutTokenT<MkFitHitWrapper> wrapperPutToken_;
0056 const edm::EDPutTokenT<MkFitClusterIndexToHit> clusterIndexPutToken_;
0057 const edm::EDPutTokenT<std::vector<float>> clusterChargePutToken_;
0058 const ConvertHitTraits convertTraits_;
0059 };
0060
0061 MkFitSiStripHitConverter::MkFitSiStripHitConverter(edm::ParameterSet const& iConfig)
0062 : stripRphiRecHitToken_{consumes<SiStripRecHit2DCollection>(iConfig.getParameter<edm::InputTag>("rphiHits"))},
0063 stripStereoRecHitToken_{consumes<SiStripRecHit2DCollection>(iConfig.getParameter<edm::InputTag>("stereoHits"))},
0064 ttrhBuilderToken_{esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(
0065 iConfig.getParameter<edm::ESInputTag>("ttrhBuilder"))},
0066 ttopoToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
0067 mkFitGeomToken_{esConsumes<MkFitGeometry, TrackerRecoGeometryRecord>()},
0068 wrapperPutToken_{produces<MkFitHitWrapper>()},
0069 clusterIndexPutToken_{produces<MkFitClusterIndexToHit>()},
0070 clusterChargePutToken_{produces<std::vector<float>>()},
0071 convertTraits_{static_cast<float>(
0072 iConfig.getParameter<edm::ParameterSet>("minGoodStripCharge").getParameter<double>("value"))} {}
0073
0074 void MkFitSiStripHitConverter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0075 edm::ParameterSetDescription desc;
0076
0077 desc.add("rphiHits", edm::InputTag{"siStripMatchedRecHits", "rphiRecHit"});
0078 desc.add("stereoHits", edm::InputTag{"siStripMatchedRecHits", "stereoRecHit"});
0079 desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
0080
0081 edm::ParameterSetDescription descCCC;
0082 descCCC.add<double>("value");
0083 desc.add("minGoodStripCharge", descCCC);
0084
0085 descriptions.add("mkFitSiStripHitConverterDefault", desc);
0086 }
0087
0088 void MkFitSiStripHitConverter::produce(edm::StreamID iID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0089 const auto& ttrhBuilder = iSetup.getData(ttrhBuilderToken_);
0090 const auto& ttopo = iSetup.getData(ttopoToken_);
0091 const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
0092
0093 MkFitHitWrapper hitWrapper;
0094 MkFitClusterIndexToHit clusterIndexToHit;
0095 std::vector<float> clusterCharge;
0096
0097 auto convert = [&](auto& hits) {
0098 return mkfit::convertHits(
0099 convertTraits_, hits, hitWrapper.hits(), clusterIndexToHit.hits(), clusterCharge, ttopo, ttrhBuilder, mkFitGeom);
0100 };
0101
0102 edm::ProductID stripClusterID;
0103 const auto& stripRphiHits = iEvent.get(stripRphiRecHitToken_);
0104 const auto& stripStereoHits = iEvent.get(stripStereoRecHitToken_);
0105 if (not stripRphiHits.empty()) {
0106 stripClusterID = convert(stripRphiHits);
0107 }
0108 if (not stripStereoHits.empty()) {
0109 auto stripStereoClusterID = convert(stripStereoHits);
0110 if (stripRphiHits.empty()) {
0111 stripClusterID = stripStereoClusterID;
0112 } else if (stripClusterID != stripStereoClusterID) {
0113 throw cms::Exception("LogicError") << "Encountered different cluster ProductIDs for strip RPhi hits ("
0114 << stripClusterID << ") and stereo (" << stripStereoClusterID << ")";
0115 }
0116 }
0117
0118 hitWrapper.setClustersID(stripClusterID);
0119
0120 iEvent.emplace(wrapperPutToken_, std::move(hitWrapper));
0121 iEvent.emplace(clusterIndexPutToken_, std::move(clusterIndexToHit));
0122 iEvent.emplace(clusterChargePutToken_, std::move(clusterCharge));
0123 }
0124
0125 DEFINE_FWK_MODULE(MkFitSiStripHitConverter);