CMS 3D CMS Logo

MkFitSiStripHitConverter.cc
Go to the documentation of this file.
2 
6 
9 
11 
13 
18 
19 #include "convertHits.h"
20 
21 namespace {
22  class ConvertHitTraits {
23  public:
24  ConvertHitTraits(float minCharge) : minGoodStripCharge_(minCharge) {}
25 
26  static constexpr bool applyCCC() { return true; }
27  static float clusterCharge(const SiStripRecHit2D& hit, DetId hitId) {
28  return siStripClusterTools::chargePerCM(hitId, hit.firstClusterRef().stripCluster());
29  }
30  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  }
34 
35  private:
36  const float minGoodStripCharge_;
37  };
38 } // namespace
39 
41 public:
42  explicit MkFitSiStripHitConverter(edm::ParameterSet const& iConfig);
43  ~MkFitSiStripHitConverter() override = default;
44 
45  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
46 
47 private:
48  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
49 
58  const ConvertHitTraits convertTraits_;
59 };
60 
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>>()},
71  convertTraits_{static_cast<float>(
72  iConfig.getParameter<edm::ParameterSet>("minGoodStripCharge").getParameter<double>("value"))} {}
73 
76 
77  desc.add("rphiHits", edm::InputTag{"siStripMatchedRecHits", "rphiRecHit"});
78  desc.add("stereoHits", edm::InputTag{"siStripMatchedRecHits", "stereoRecHit"});
79  desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
80 
82  descCCC.add<double>("value");
83  desc.add("minGoodStripCharge", descCCC);
84 
85  descriptions.add("mkFitSiStripHitConverterDefault", desc);
86 }
87 
89  const auto& ttrhBuilder = iSetup.getData(ttrhBuilderToken_);
90  const auto& ttopo = iSetup.getData(ttopoToken_);
91  const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
92 
93  MkFitHitWrapper hitWrapper;
94  MkFitClusterIndexToHit clusterIndexToHit;
95  std::vector<float> clusterCharge;
96 
97  auto convert = [&](auto& hits) {
98  return mkfit::convertHits(
99  convertTraits_, hits, hitWrapper.hits(), clusterIndexToHit.hits(), clusterCharge, ttopo, ttrhBuilder, mkFitGeom);
100  };
101 
102  edm::ProductID stripClusterID;
103  const auto& stripRphiHits = iEvent.get(stripRphiRecHitToken_);
104  const auto& stripStereoHits = iEvent.get(stripStereoRecHitToken_);
105  if (not stripRphiHits.empty()) {
106  stripClusterID = convert(stripRphiHits);
107  }
108  if (not stripStereoHits.empty()) {
109  auto stripStereoClusterID = convert(stripStereoHits);
110  if (stripRphiHits.empty()) {
111  stripClusterID = stripStereoClusterID;
112  } else if (stripClusterID != stripStereoClusterID) {
113  throw cms::Exception("LogicError") << "Encountered different cluster ProductIDs for strip RPhi hits ("
114  << stripClusterID << ") and stereo (" << stripStereoClusterID << ")";
115  }
116  }
117 
118  hitWrapper.setClustersID(stripClusterID);
119 
120  iEvent.emplace(wrapperPutToken_, std::move(hitWrapper));
121  iEvent.emplace(clusterIndexPutToken_, std::move(clusterIndexToHit));
122  iEvent.emplace(clusterChargePutToken_, std::move(clusterCharge));
123 }
124 
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::EDGetTokenT< SiStripRecHit2DCollection > stripStereoRecHitToken_
void setupAsStrip(unsigned int id, int cpcm, int rows)
Definition: Hit.h:244
const edm::ESGetToken< MkFitGeometry, TrackerRecoGeometryRecord > mkFitGeomToken_
const edm::EDGetTokenT< SiStripRecHit2DCollection > stripRphiRecHitToken_
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
float chargePerCM(DetId detid, Iter a, Iter b)
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
edm::ProductID convertHits(const Traits &traits, const HitCollection &hits, mkfit::HitVec &mkFitHits, std::vector< TrackingRecHit const *> &clusterIndexToHit, std::vector< float > &clusterChargeVec, const TrackerTopology &ttopo, const TransientTrackingRecHitBuilder &ttrhBuilder, const MkFitGeometry &mkFitGeom)
Definition: convertHits.h:25
const edm::EDPutTokenT< MkFitClusterIndexToHit > clusterIndexPutToken_
int iEvent
Definition: GenABIO.cc:224
const edm::EDPutTokenT< std::vector< float > > clusterChargePutToken_
SiStripCluster const & amplitudes() const
auto size() const
const edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > ttrhBuilderToken_
const edm::EDPutTokenT< MkFitHitWrapper > wrapperPutToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MkFitSiStripHitConverter(edm::ParameterSet const &iConfig)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const ConvertHitTraits convertTraits_
W convert(V v)
Definition: ExtVec.h:66
Definition: DetId.h:17
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void setClustersID(edm::ProductID id)
mkfit::HitVec & hits()
~MkFitSiStripHitConverter() override=default
std::vector< TrackingRecHit const * > & hits()
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:511