CMS 3D CMS Logo

PixelVertexProducerCUDA.cc
Go to the documentation of this file.
1 #include <cuda_runtime.h>
2 
19 
20 #include "gpuVertexFinder.h"
21 
22 #undef PIXVERTEX_DEBUG_PRODUCE
23 
25 public:
26  explicit PixelVertexProducerCUDA(const edm::ParameterSet& iConfig);
27  ~PixelVertexProducerCUDA() override = default;
28 
29  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
30 
31 private:
32  void produceOnGPU(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const;
33  void produceOnCPU(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const;
34  void produce(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
35 
36  bool onGPU_;
37 
42 
44 
45  // Tracking cuts before sending tracks to vertex algo
46  const float ptMin_;
47  const float ptMax_;
48 };
49 
51  : onGPU_(conf.getParameter<bool>("onGPU")),
52  gpuAlgo_(conf.getParameter<bool>("oneKernel"),
53  conf.getParameter<bool>("useDensity"),
54  conf.getParameter<bool>("useDBSCAN"),
55  conf.getParameter<bool>("useIterative"),
56  conf.getParameter<int>("minT"),
57  conf.getParameter<double>("eps"),
58  conf.getParameter<double>("errmax"),
59  conf.getParameter<double>("chi2max")),
60  ptMin_(conf.getParameter<double>("PtMin")), // 0.5 GeV
61  ptMax_(conf.getParameter<double>("PtMax")) // 75. GeV
62 {
63  if (onGPU_) {
65  consumes<cms::cuda::Product<PixelTrackHeterogeneous>>(conf.getParameter<edm::InputTag>("pixelTrackSrc"));
66  tokenGPUVertex_ = produces<ZVertexCUDAProduct>();
67  } else {
68  tokenCPUTrack_ = consumes<PixelTrackHeterogeneous>(conf.getParameter<edm::InputTag>("pixelTrackSrc"));
69  tokenCPUVertex_ = produces<ZVertexHeterogeneous>();
70  }
71 }
72 
75 
76  // Only one of these three algos can be used at once.
77  // Maybe this should become a Plugin Factory
78  desc.add<bool>("onGPU", true);
79  desc.add<bool>("oneKernel", true);
80  desc.add<bool>("useDensity", true);
81  desc.add<bool>("useDBSCAN", false);
82  desc.add<bool>("useIterative", false);
83 
84  desc.add<int>("minT", 2); // min number of neighbours to be "core"
85  desc.add<double>("eps", 0.07); // max absolute distance to cluster
86  desc.add<double>("errmax", 0.01); // max error to be "seed"
87  desc.add<double>("chi2max", 9.); // max normalized distance to cluster
88 
89  desc.add<double>("PtMin", 0.5);
90  desc.add<double>("PtMax", 75.);
91  desc.add<edm::InputTag>("pixelTrackSrc", edm::InputTag("pixelTracksCUDA"));
92 
93  auto label = "pixelVerticesCUDA";
94  descriptions.add(label, desc);
95 }
96 
99  const edm::EventSetup& iSetup) const {
101  iEvent.getByToken(tokenGPUTrack_, hTracks);
102 
103  cms::cuda::ScopedContextProduce ctx{*hTracks};
104  auto const* tracks = ctx.get(*hTracks).get();
105 
106  assert(tracks);
107 
108  ctx.emplace(iEvent, tokenGPUVertex_, gpuAlgo_.makeAsync(ctx.stream(), tracks, ptMin_, ptMax_));
109 }
110 
113  const edm::EventSetup& iSetup) const {
114  auto const* tracks = iEvent.get(tokenCPUTrack_).get();
115  assert(tracks);
116 
117 #ifdef PIXVERTEX_DEBUG_PRODUCE
118  auto const& tsoa = *tracks;
119  auto maxTracks = tsoa.stride();
120  std::cout << "size of SoA " << sizeof(tsoa) << " stride " << maxTracks << std::endl;
121 
122  int32_t nt = 0;
123  for (int32_t it = 0; it < maxTracks; ++it) {
124  auto nHits = tsoa.nHits(it);
125  assert(nHits == int(tsoa.hitIndices.size(it)));
126  if (nHits == 0)
127  break; // this is a guard: maybe we need to move to nTracks...
128  nt++;
129  }
130  std::cout << "found " << nt << " tracks in cpu SoA for Vertexing at " << tracks << std::endl;
131 #endif // PIXVERTEX_DEBUG_PRODUCE
132 
134 }
135 
137  if (onGPU_) {
138  produceOnGPU(streamID, iEvent, iSetup);
139  } else {
140  produceOnCPU(streamID, iEvent, iSetup);
141  }
142 }
143 
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
ZVertexHeterogeneous makeAsync(cudaStream_t stream, TkSoA const *tksoa, float ptMin, float ptMax) const
PixelVertexProducerCUDA(const edm::ParameterSet &iConfig)
edm::EDPutTokenT< ZVertexCUDAProduct > tokenGPUVertex_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDPutTokenT< ZVertexHeterogeneous > tokenCPUVertex_
edm::EDGetTokenT< cms::cuda::Product< PixelTrackHeterogeneous > > tokenGPUTrack_
assert(be >=bs)
char const * label
void produceOnGPU(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const
int iEvent
Definition: GenABIO.cc:224
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produceOnCPU(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const
void produce(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
int nt
Definition: AMPTWrapper.h:42
auto const & tracks
cannot be loose
void add(std::string const &label, ParameterSetDescription const &psetDescription)
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
~PixelVertexProducerCUDA() override=default
ZVertexHeterogeneous make(TkSoA const *tksoa, float ptMin, float ptMax) const
edm::EDGetTokenT< PixelTrackHeterogeneous > tokenCPUTrack_
const gpuVertexFinder::Producer gpuAlgo_