CMS 3D CMS Logo

MkFitSiPixelHitConverter.cc
Go to the documentation of this file.
2 
6 
8 
10 
12 
17 
18 #include "convertHits.h"
19 
20 namespace {
21  struct ConvertHitTraits {
22  static constexpr bool applyCCC() { return false; }
23  static std::nullptr_t clusterCharge(const SiPixelRecHit& hit, DetId hitId) { return nullptr; }
24  static bool passCCC(std::nullptr_t) { return true; }
25  static void setDetails(mkfit::Hit& mhit, const SiPixelCluster& cluster, int shortId, std::nullptr_t) {
26  mhit.setupAsPixel(shortId, cluster.sizeX(), cluster.sizeY());
27  }
28  };
29 } // namespace
30 
32 public:
33  explicit MkFitSiPixelHitConverter(edm::ParameterSet const& iConfig);
34  ~MkFitSiPixelHitConverter() override = default;
35 
36  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
37 
38 private:
39  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
40 
47 };
48 
50  : pixelRecHitToken_{consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("hits"))},
51  ttrhBuilderToken_{esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(
52  iConfig.getParameter<edm::ESInputTag>("ttrhBuilder"))},
53  ttopoToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
54  mkFitGeomToken_{esConsumes<MkFitGeometry, TrackerRecoGeometryRecord>()},
55  wrapperPutToken_{produces<MkFitHitWrapper>()},
56  clusterIndexPutToken_{produces<MkFitClusterIndexToHit>()} {}
57 
60 
61  desc.add("hits", edm::InputTag{"siPixelRecHits"});
62  desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
63 
64  descriptions.addWithDefaultLabel(desc);
65 }
66 
68  const auto& ttrhBuilder = iSetup.getData(ttrhBuilderToken_);
69  const auto& ttopo = iSetup.getData(ttopoToken_);
70  const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
71 
72  MkFitHitWrapper hitWrapper;
73  MkFitClusterIndexToHit clusterIndexToHit;
74 
75  std::vector<float> dummy;
76  auto pixelClusterID = mkfit::convertHits(ConvertHitTraits{},
78  hitWrapper.hits(),
79  clusterIndexToHit.hits(),
80  dummy,
81  ttopo,
82  ttrhBuilder,
83  mkFitGeom);
84 
85  hitWrapper.setClustersID(pixelClusterID);
86 
87  iEvent.emplace(wrapperPutToken_, std::move(hitWrapper));
88  iEvent.emplace(clusterIndexPutToken_, std::move(clusterIndexToHit));
89 }
90 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
~MkFitSiPixelHitConverter() override=default
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int sizeY() const
const edm::EDPutTokenT< MkFitClusterIndexToHit > clusterIndexPutToken_
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
void setupAsPixel(unsigned int id, int rows, int cols)
Definition: Hit.h:237
const edm::EDPutTokenT< MkFitHitWrapper > wrapperPutToken_
int iEvent
Definition: GenABIO.cc:224
int sizeX() const
const edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > ttrhBuilderToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool getData(T &iHolder) const
Definition: EventSetup.h:122
Definition: DetId.h:17
const edm::ESGetToken< MkFitGeometry, TrackerRecoGeometryRecord > mkFitGeomToken_
void setClustersID(edm::ProductID id)
mkfit::HitVec & hits()
Pixel cluster – collection of neighboring pixels above threshold.
MkFitSiPixelHitConverter(edm::ParameterSet const &iConfig)
const edm::EDGetTokenT< SiPixelRecHitCollection > pixelRecHitToken_
std::vector< TrackingRecHit const * > & hits()
def move(src, dest)
Definition: eostools.py:511
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
Our base class.
Definition: SiPixelRecHit.h:23