CMS 3D CMS Logo

MkFitPhase2HitConverter.cc
Go to the documentation of this file.
2 
6 
8 
10 
12 
17 
18 #include "convertHits.h"
19 
20 namespace {
21  class ConvertHitTraitsPhase2 {
22  public:
23  static constexpr bool applyCCC() { return false; }
24  static std::nullptr_t clusterCharge(const Phase2TrackerRecHit1D& hit, DetId hitId) { return nullptr; }
25  static bool passCCC(std::nullptr_t) { return true; }
26  static void setDetails(mkfit::Hit& mhit, const Phase2TrackerCluster1D& cluster, int shortId, std::nullptr_t) {
27  mhit.setupAsStrip(shortId, (1 << 8) - 1, cluster.size());
28  }
29  };
30 } // namespace
31 
33 public:
34  explicit MkFitPhase2HitConverter(edm::ParameterSet const& iConfig);
35  ~MkFitPhase2HitConverter() override = default;
36 
37  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
38 
39 private:
40  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
41 
49  const ConvertHitTraitsPhase2 convertTraits_;
50 };
51 
53  : siPhase2RecHitToken_{consumes<Phase2TrackerRecHit1DCollectionNew>(
54  iConfig.getParameter<edm::InputTag>("siPhase2Hits"))},
55  ttrhBuilderToken_{esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(
56  iConfig.getParameter<edm::ESInputTag>("ttrhBuilder"))},
57  ttopoToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
58  mkFitGeomToken_{esConsumes<MkFitGeometry, TrackerRecoGeometryRecord>()},
59  wrapperPutToken_{produces<MkFitHitWrapper>()},
60  clusterIndexPutToken_{produces<MkFitClusterIndexToHit>()},
61  clusterChargePutToken_{produces<std::vector<float>>()},
62  convertTraits_{} {}
63 
66 
67  desc.add("siPhase2Hits", edm::InputTag{"siPhase2RecHits"});
68  desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
69 
70  descriptions.add("mkFitPhase2HitConverterDefault", desc);
71 }
72 
74  const auto& ttrhBuilder = iSetup.getData(ttrhBuilderToken_);
75  const auto& ttopo = iSetup.getData(ttopoToken_);
76  const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
77 
78  MkFitHitWrapper hitWrapper;
79  MkFitClusterIndexToHit clusterIndexToHit;
80  std::vector<float> clusterCharge;
81 
82  edm::ProductID stripClusterID;
83  const auto& phase2Hits = iEvent.get(siPhase2RecHitToken_);
84  std::vector<float> dummy;
85  if (not phase2Hits.empty()) {
86  stripClusterID = mkfit::convertHits(ConvertHitTraitsPhase2{},
87  phase2Hits,
88  hitWrapper.hits(),
89  clusterIndexToHit.hits(),
90  dummy,
91  ttopo,
92  ttrhBuilder,
93  mkFitGeom);
94  }
95 
96  hitWrapper.setClustersID(stripClusterID);
97 
98  iEvent.emplace(wrapperPutToken_, std::move(hitWrapper));
99  iEvent.emplace(clusterIndexPutToken_, std::move(clusterIndexToHit));
100  iEvent.emplace(clusterChargePutToken_, std::move(clusterCharge));
101 }
102 
const edm::EDPutTokenT< std::vector< float > > clusterChargePutToken_
const edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > ttrhBuilderToken_
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
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
void setupAsStrip(unsigned int id, int cpcm, int rows)
Definition: Hit.h:244
const ConvertHitTraitsPhase2 convertTraits_
const edm::EDGetTokenT< Phase2TrackerRecHit1DCollectionNew > siPhase2RecHitToken_
~MkFitPhase2HitConverter() override=default
const edm::EDPutTokenT< MkFitHitWrapper > wrapperPutToken_
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::ESGetToken< MkFitGeometry, TrackerRecoGeometryRecord > mkFitGeomToken_
int iEvent
Definition: GenABIO.cc:224
const edm::EDPutTokenT< MkFitClusterIndexToHit > clusterIndexPutToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Definition: DetId.h:17
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void setClustersID(edm::ProductID id)
mkfit::HitVec & hits()
MkFitPhase2HitConverter(edm::ParameterSet const &iConfig)
std::vector< TrackingRecHit const * > & hits()
def move(src, dest)
Definition: eostools.py:511