@@ -22,14 +22,16 @@ namespace {
2222 class ConvertHitTraits {
2323 public:
2424 ConvertHitTraits (float minCharge) : minGoodStripCharge_(minCharge) {}
25+ using Clusters = edmNew::DetSetVector<SiStripCluster>;
26+ using Cluster = Clusters::data_type;
2527
2628 static constexpr bool applyCCC () { return true ; }
27- static float clusterCharge ( const SiStripRecHit2D& hit, DetId hitId) {
28- return siStripClusterTools::chargePerCM (hitId, hit. firstClusterRef (). stripCluster ());
29- }
29+ static float chargeScale (DetId id) { return siStripClusterTools::sensorThicknessInverse (id); }
30+ static const Cluster& cluster ( const Clusters& prod, unsigned int index) { return prod. data ()[index]; }
31+ static float clusterCharge ( const Cluster& clu, float scale) { return clu. charge () * scale; }
3032 bool passCCC (float charge) const { return charge > minGoodStripCharge_; }
31- static void setDetails (mkfit::Hit& mhit, const SiStripCluster& cluster , int shortId, float charge) {
32- mhit.setupAsStrip (shortId, charge, cluster. amplitudes () .size ());
33+ static void setDetails (mkfit::Hit& mhit, const Cluster& clu , int shortId, float charge) {
34+ mhit.setupAsStrip (shortId, charge, clu .size ());
3335 }
3436
3537 private:
@@ -49,25 +51,28 @@ class MkFitSiStripHitConverter : public edm::global::EDProducer<> {
4951
5052 const edm::EDGetTokenT<SiStripRecHit2DCollection> stripRphiRecHitToken_;
5153 const edm::EDGetTokenT<SiStripRecHit2DCollection> stripStereoRecHitToken_;
54+ const edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster>> stripClusterToken_;
5255 const edm::ESGetToken<TransientTrackingRecHitBuilder, TransientRecHitRecord> ttrhBuilderToken_;
5356 const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> ttopoToken_;
5457 const edm::ESGetToken<MkFitGeometry, TrackerRecoGeometryRecord> mkFitGeomToken_;
5558 const edm::EDPutTokenT<MkFitHitWrapper> wrapperPutToken_;
5659 const edm::EDPutTokenT<MkFitClusterIndexToHit> clusterIndexPutToken_;
60+ const edm::EDPutTokenT<std::vector<int >> layerIndexPutToken_;
5761 const edm::EDPutTokenT<std::vector<float >> clusterChargePutToken_;
5862 const ConvertHitTraits convertTraits_;
5963};
6064
6165MkFitSiStripHitConverter::MkFitSiStripHitConverter (edm::ParameterSet const & iConfig)
62- : stripRphiRecHitToken_{consumes<SiStripRecHit2DCollection>(iConfig.getParameter <edm::InputTag>(" rphiHits" ))},
63- stripStereoRecHitToken_{consumes<SiStripRecHit2DCollection>(iConfig.getParameter <edm::InputTag>(" stereoHits" ))},
64- ttrhBuilderToken_{esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(
65- iConfig.getParameter <edm::ESInputTag>(" ttrhBuilder" ))},
66- ttopoToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
67- mkFitGeomToken_{esConsumes<MkFitGeometry, TrackerRecoGeometryRecord>()},
68- wrapperPutToken_{produces<MkFitHitWrapper>()},
69- clusterIndexPutToken_{produces<MkFitClusterIndexToHit>()},
70- clusterChargePutToken_{produces<std::vector<float >>()},
66+ : stripRphiRecHitToken_{consumes (iConfig.getParameter <edm::InputTag>(" rphiHits" ))},
67+ stripStereoRecHitToken_{consumes (iConfig.getParameter <edm::InputTag>(" stereoHits" ))},
68+ stripClusterToken_{consumes (iConfig.getParameter <edm::InputTag>(" clusters" ))},
69+ ttrhBuilderToken_{esConsumes (iConfig.getParameter <edm::ESInputTag>(" ttrhBuilder" ))},
70+ ttopoToken_{esConsumes ()},
71+ mkFitGeomToken_{esConsumes ()},
72+ wrapperPutToken_{produces ()},
73+ clusterIndexPutToken_{produces ()},
74+ layerIndexPutToken_{produces ()},
75+ clusterChargePutToken_{produces ()},
7176 convertTraits_{static_cast <float >(
7277 iConfig.getParameter <edm::ParameterSet>(" minGoodStripCharge" ).getParameter <double >(" value" ))} {}
7378
@@ -76,6 +81,7 @@ void MkFitSiStripHitConverter::fillDescriptions(edm::ConfigurationDescriptions&
7681
7782 desc.add (" rphiHits" , edm::InputTag{" siStripMatchedRecHits" , " rphiRecHit" });
7883 desc.add (" stereoHits" , edm::InputTag{" siStripMatchedRecHits" , " stereoRecHit" });
84+ desc.add (" clusters" , edm::InputTag{" siStripClusters" });
7985 desc.add (" ttrhBuilder" , edm::ESInputTag{" " , " WithTrackAngle" });
8086
8187 edm::ParameterSetDescription descCCC;
@@ -92,16 +98,29 @@ void MkFitSiStripHitConverter::produce(edm::StreamID iID, edm::Event& iEvent, co
9298
9399 MkFitHitWrapper hitWrapper;
94100 MkFitClusterIndexToHit clusterIndexToHit;
101+ std::vector<int > layerIndexToHit;
95102 std::vector<float > clusterCharge;
96103
97- auto convert = [&](auto & hits) {
98- return mkfit::convertHits (
99- convertTraits_, hits, hitWrapper.hits (), clusterIndexToHit.hits (), clusterCharge, ttopo, ttrhBuilder, mkFitGeom);
100- };
101-
102104 edm::ProductID stripClusterID;
103105 const auto & stripRphiHits = iEvent.get (stripRphiRecHitToken_);
104106 const auto & stripStereoHits = iEvent.get (stripStereoRecHitToken_);
107+ const auto maxSizeGuess (stripRphiHits.dataSize () + stripStereoHits.dataSize ());
108+ auto const & clusters = iEvent.get (stripClusterToken_);
109+
110+ auto convert = [&](auto & hits) {
111+ return mkfit::convertHits (convertTraits_,
112+ hits,
113+ clusters,
114+ hitWrapper.hits (),
115+ clusterIndexToHit.hits (),
116+ layerIndexToHit,
117+ clusterCharge,
118+ ttopo,
119+ ttrhBuilder,
120+ mkFitGeom,
121+ maxSizeGuess);
122+ };
123+
105124 if (not stripRphiHits.empty ()) {
106125 stripClusterID = convert (stripRphiHits);
107126 }
@@ -119,6 +138,7 @@ void MkFitSiStripHitConverter::produce(edm::StreamID iID, edm::Event& iEvent, co
119138
120139 iEvent.emplace (wrapperPutToken_, std::move (hitWrapper));
121140 iEvent.emplace (clusterIndexPutToken_, std::move (clusterIndexToHit));
141+ iEvent.emplace (layerIndexPutToken_, std::move (layerIndexToHit));
122142 iEvent.emplace (clusterChargePutToken_, std::move (clusterCharge));
123143}
124144
0 commit comments