CMS 3D CMS Logo

PixelVertexProducerCUDA.cc
Go to the documentation of this file.
1 #include <cuda_runtime.h>
2 
20 
25 
26 #include "gpuVertexFinder.h"
27 
28 #undef PIXVERTEX_DEBUG_PRODUCE
29 
30 template <typename TrackerTraits>
35 
36 public:
37  explicit PixelVertexProducerCUDAT(const edm::ParameterSet& iConfig);
38  ~PixelVertexProducerCUDAT() override = default;
39 
40  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
41 
42 private:
43  void produceOnGPU(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const;
44  void produceOnCPU(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const;
45  void produce(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
46 
47  bool onGPU_;
48 
53 
55 
56  // Tracking cuts before sending tracks to vertex algo
57  const float ptMin_;
58  const float ptMax_;
59 };
60 
61 template <typename TrackerTraits>
63  : onGPU_(conf.getParameter<bool>("onGPU")),
64  gpuAlgo_(conf.getParameter<bool>("oneKernel"),
65  conf.getParameter<bool>("useDensity"),
66  conf.getParameter<bool>("useDBSCAN"),
67  conf.getParameter<bool>("useIterative"),
68  conf.getParameter<int>("minT"),
69  conf.getParameter<double>("eps"),
70  conf.getParameter<double>("errmax"),
71  conf.getParameter<double>("chi2max")),
72  ptMin_(conf.getParameter<double>("PtMin")), // 0.5 GeV
73  ptMax_(conf.getParameter<double>("PtMax")) // 75. GeV
74 {
75  if (onGPU_) {
76  tokenGPUTrack_ = consumes(conf.getParameter<edm::InputTag>("pixelTrackSrc"));
78  } else {
79  tokenCPUTrack_ = consumes(conf.getParameter<edm::InputTag>("pixelTrackSrc"));
81  }
82 }
83 
84 template <typename TrackerTraits>
87 
88  // Only one of these three algos can be used at once.
89  // Maybe this should become a Plugin Factory
90  desc.add<bool>("onGPU", true);
91  desc.add<bool>("oneKernel", true);
92  desc.add<bool>("useDensity", true);
93  desc.add<bool>("useDBSCAN", false);
94  desc.add<bool>("useIterative", false);
95 
96  desc.add<int>("minT", 2); // min number of neighbours to be "core"
97  desc.add<double>("eps", 0.07); // max absolute distance to cluster
98  desc.add<double>("errmax", 0.01); // max error to be "seed"
99  desc.add<double>("chi2max", 9.); // max normalized distance to cluster
100 
101  desc.add<double>("PtMin", 0.5);
102  desc.add<double>("PtMax", 75.);
103  desc.add<edm::InputTag>("pixelTrackSrc", edm::InputTag("pixelTracksCUDA"));
104 
105  descriptions.addWithDefaultLabel(desc);
106 }
107 
108 template <typename TrackerTraits>
111  const edm::EventSetup& iSetup) const {
113  auto hTracks = iEvent.getHandle(tokenGPUTrack_);
114 
115  cms::cuda::ScopedContextProduce ctx{*hTracks};
116  auto& tracks = ctx.get(*hTracks);
117 
118  ctx.emplace(iEvent, tokenGPUVertex_, gpuAlgo_.makeAsync(ctx.stream(), tracks.view(), ptMin_, ptMax_));
119 }
120 
121 template <typename TrackerTraits>
124  const edm::EventSetup& iSetup) const {
125  auto& tracks = iEvent.get(tokenCPUTrack_);
126 
127 #ifdef PIXVERTEX_DEBUG_PRODUCE
128  auto const& tsoa = *tracks;
129  auto maxTracks = tsoa.stride();
130  std::cout << "size of SoA " << sizeof(tsoa) << " stride " << maxTracks << std::endl;
131 
132  int32_t nt = 0;
133  for (int32_t it = 0; it < maxTracks; ++it) {
135  assert(nHits == int(tracks.view().hitIndices().size(it)));
136  if (nHits == 0)
137  break; // this is a guard: maybe we need to move to nTracks...
138  nt++;
139  }
140  std::cout << "found " << nt << " tracks in cpu SoA for Vertexing at " << tracks << std::endl;
141 #endif // PIXVERTEX_DEBUG_PRODUCE
142 
143  iEvent.emplace(tokenCPUVertex_, gpuAlgo_.make(tracks.view(), ptMin_, ptMax_));
144 }
145 
146 template <typename TrackerTraits>
149  const edm::EventSetup& iSetup) const {
150  if (onGPU_) {
151  produceOnGPU(streamID, iEvent, iSetup);
152  } else {
153  produceOnCPU(streamID, iEvent, iSetup);
154  }
155 }
156 
159 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< TracksSoAHost > tokenCPUTrack_
void produceOnGPU(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const
~PixelVertexProducerCUDAT() override=default
assert(be >=bs)
edm::EDPutTokenT< cms::cuda::Product< ZVertexSoADevice > > tokenGPUVertex_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
edm::EDPutTokenT< ZVertexSoAHost > tokenCPUVertex_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produceOnCPU(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const
edm::EDGetTokenT< cms::cuda::Product< TracksSoADevice > > tokenGPUTrack_
int nt
Definition: AMPTWrapper.h:42
static constexpr __host__ __device__ int nHits(const TrackSoAConstView &tracks, int i)
maxTracks
Definition: DMR_cfg.py:158
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
PixelVertexProducerCUDAT(const edm::ParameterSet &iConfig)
void produce(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits