CMS 3D CMS Logo

CAHitNtupletGeneratorOnGPU.cc
Go to the documentation of this file.
1 //
2 // Original Author: Felice Pantaleo, CERN
3 //
4 
5 // #define GPU_DEBUG
6 // #define DUMP_GPU_TK_TUPLES
7 
8 #include <array>
9 #include <cassert>
10 #include <functional>
11 #include <vector>
12 
27 
29 
30 namespace {
31 
32  using namespace caHitNtupletGenerator;
33  using namespace gpuPixelDoublets;
34  using namespace pixelTopology;
35  using namespace pixelTrack;
36 
37  template <typename T>
38  T sqr(T x) {
39  return x * x;
40  }
41 
42  //Common Params
43  AlgoParams makeCommonParams(edm::ParameterSet const& cfg) {
44  return AlgoParams({cfg.getParameter<bool>("onGPU"),
45  cfg.getParameter<unsigned int>("minHitsForSharingCut"),
46  cfg.getParameter<bool>("useRiemannFit"),
47  cfg.getParameter<bool>("fitNas4"),
48  cfg.getParameter<bool>("includeJumpingForwardDoublets"),
49  cfg.getParameter<bool>("earlyFishbone"),
50  cfg.getParameter<bool>("lateFishbone"),
51  cfg.getParameter<bool>("fillStatistics"),
52  cfg.getParameter<bool>("doSharedHitCut"),
53  cfg.getParameter<bool>("dupPassThrough"),
54  cfg.getParameter<bool>("useSimpleTripletCleaner")});
55  }
56 
57  //This is needed to have the partial specialization for isPhase1Topology/isPhase2Topology
58  template <typename TrackerTraits, typename Enable = void>
59  struct topologyCuts {};
60 
61  template <typename TrackerTraits>
62  struct topologyCuts<TrackerTraits, isPhase1Topology<TrackerTraits>> {
63  static constexpr CAParamsT<TrackerTraits> makeCACuts(edm::ParameterSet const& cfg) {
65  cfg.getParameter<unsigned int>("maxNumberOfDoublets"),
66  cfg.getParameter<unsigned int>("minHitsPerNtuplet"),
67  static_cast<float>(cfg.getParameter<double>("ptmin")),
68  static_cast<float>(cfg.getParameter<double>("CAThetaCutBarrel")),
69  static_cast<float>(cfg.getParameter<double>("CAThetaCutForward")),
70  static_cast<float>(cfg.getParameter<double>("hardCurvCut")),
71  static_cast<float>(cfg.getParameter<double>("dcaCutInnerTriplet")),
72  static_cast<float>(cfg.getParameter<double>("dcaCutOuterTriplet")),
73  }};
74  };
75 
76  static constexpr pixelTrack::QualityCutsT<TrackerTraits> makeQualityCuts(edm::ParameterSet const& pset) {
77  auto coeff = pset.getParameter<std::array<double, 2>>("chi2Coeff");
78  auto ptMax = pset.getParameter<double>("chi2MaxPt");
79 
80  coeff[1] = (coeff[1] - coeff[0]) / log2(ptMax);
81  return pixelTrack::QualityCutsT<TrackerTraits>{// polynomial coefficients for the pT-dependent chi2 cut
82  {(float)coeff[0], (float)coeff[1], 0.f, 0.f},
83  // max pT used to determine the chi2 cut
84  (float)ptMax,
85  // chi2 scale factor: 8 for broken line fit, ?? for Riemann fit
86  (float)pset.getParameter<double>("chi2Scale"),
87  // regional cuts for triplets
88  {(float)pset.getParameter<double>("tripletMaxTip"),
89  (float)pset.getParameter<double>("tripletMinPt"),
90  (float)pset.getParameter<double>("tripletMaxZip")},
91  // regional cuts for quadruplets
92  {(float)pset.getParameter<double>("quadrupletMaxTip"),
93  (float)pset.getParameter<double>("quadrupletMinPt"),
94  (float)pset.getParameter<double>("quadrupletMaxZip")}};
95  }
96  };
97 
98  template <typename TrackerTraits>
99  struct topologyCuts<TrackerTraits, isPhase2Topology<TrackerTraits>> {
100  static constexpr CAParamsT<TrackerTraits> makeCACuts(edm::ParameterSet const& cfg) {
101  return CAParamsT<TrackerTraits>{{cfg.getParameter<unsigned int>("maxNumberOfDoublets"),
102  cfg.getParameter<unsigned int>("minHitsPerNtuplet"),
103  static_cast<float>(cfg.getParameter<double>("ptmin")),
104  static_cast<float>(cfg.getParameter<double>("CAThetaCutBarrel")),
105  static_cast<float>(cfg.getParameter<double>("CAThetaCutForward")),
106  static_cast<float>(cfg.getParameter<double>("hardCurvCut")),
107  static_cast<float>(cfg.getParameter<double>("dcaCutInnerTriplet")),
108  static_cast<float>(cfg.getParameter<double>("dcaCutOuterTriplet"))},
109  {(bool)cfg.getParameter<bool>("includeFarForwards")}};
110  }
111 
112  static constexpr pixelTrack::QualityCutsT<TrackerTraits> makeQualityCuts(edm::ParameterSet const& pset) {
114  static_cast<float>(pset.getParameter<double>("maxChi2")),
115  static_cast<float>(pset.getParameter<double>("minPt")),
116  static_cast<float>(pset.getParameter<double>("maxTip")),
117  static_cast<float>(pset.getParameter<double>("maxZip")),
118  };
119  }
120  };
121 
122  //Cell Cuts, as they are the cuts have the same logic for Phase2 and Phase1
123  //keeping them separate would allow further differentiation in the future
124  //moving them to topologyCuts and using the same syntax
125  template <typename TrackerTraits>
126  CellCutsT<TrackerTraits> makeCellCuts(edm::ParameterSet const& cfg) {
127  return CellCutsT<TrackerTraits>{cfg.getParameter<bool>("doClusterCut"),
128  cfg.getParameter<bool>("doZ0Cut"),
129  cfg.getParameter<bool>("doPtCut"),
130  cfg.getParameter<bool>("idealConditions"),
131  (float)cfg.getParameter<double>("z0Cut"),
132  (float)cfg.getParameter<double>("ptCut"),
133  cfg.getParameter<std::vector<int>>("phiCuts")};
134  }
135 
136 } // namespace
137 
138 using namespace std;
139 
140 template <typename TrackerTraits>
143  : m_params(makeCommonParams(cfg),
144  makeCellCuts<TrackerTraits>(cfg),
145  topologyCuts<TrackerTraits>::makeQualityCuts(cfg.getParameterSet("trackQualityCuts")),
146  topologyCuts<TrackerTraits>::makeCACuts(cfg)) {
147 #ifdef DUMP_GPU_TK_TUPLES
148  printf("TK: %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s\n",
149  "tid",
150  "qual",
151  "nh",
152  "nl",
153  "charge",
154  "pt",
155  "eta",
156  "phi",
157  "tip",
158  "zip",
159  "chi2",
160  "h1",
161  "h2",
162  "h3",
163  "h4",
164  "h5",
165  "hn");
166 #endif
167 }
168 
169 template <typename TrackerTraits>
171  fillDescriptionsCommon(desc);
172  edm::LogWarning("CAHitNtupletGeneratorOnGPU::fillDescriptions")
173  << "Note: this fillDescriptions is a dummy one. Most probably you are missing some parameters. \n"
174  "please implement your TrackerTraits descriptions in CAHitNtupletGeneratorOnGPU. \n";
175 }
176 
177 template <>
179  fillDescriptionsCommon(desc);
180 
181  desc.add<bool>("idealConditions", true);
182  desc.add<bool>("includeJumpingForwardDoublets", false);
183  desc.add<double>("z0Cut", 12.0);
184  desc.add<double>("ptCut", 0.5);
185 
187  trackQualityCuts.add<double>("chi2MaxPt", 10.)->setComment("max pT used to determine the pT-dependent chi2 cut");
188  trackQualityCuts.add<std::vector<double>>("chi2Coeff", {0.9, 1.8})->setComment("chi2 at 1GeV and at ptMax above");
189  trackQualityCuts.add<double>("chi2Scale", 8.)
190  ->setComment(
191  "Factor to multiply the pT-dependent chi2 cut (currently: 8 for the broken line fit, ?? for the Riemann "
192  "fit)");
193  trackQualityCuts.add<double>("tripletMinPt", 0.5)->setComment("Min pT for triplets, in GeV");
194  trackQualityCuts.add<double>("tripletMaxTip", 0.3)->setComment("Max |Tip| for triplets, in cm");
195  trackQualityCuts.add<double>("tripletMaxZip", 12.)->setComment("Max |Zip| for triplets, in cm");
196  trackQualityCuts.add<double>("quadrupletMinPt", 0.3)->setComment("Min pT for quadruplets, in GeV");
197  trackQualityCuts.add<double>("quadrupletMaxTip", 0.5)->setComment("Max |Tip| for quadruplets, in cm");
198  trackQualityCuts.add<double>("quadrupletMaxZip", 12.)->setComment("Max |Zip| for quadruplets, in cm");
199 
200  desc.add<std::vector<int>>(
201  "phiCuts", std::vector<int>(std::begin(phase1PixelTopology::phicuts), std::end(phase1PixelTopology::phicuts)))
202  ->setComment("Cuts in phi for cells");
203 
204  desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
205  ->setComment(
206  "Quality cuts based on the results of the track fit:\n - apply a pT-dependent chi2 cut;\n - apply \"region "
207  "cuts\" based on the fit results (pT, Tip, Zip).");
208 }
209 
210 template <>
212  fillDescriptionsCommon(desc);
213 
214  desc.add<bool>("idealConditions", false);
215  desc.add<bool>("includeJumpingForwardDoublets", false);
216  desc.add<double>("z0Cut", 10.0);
217  desc.add<double>("ptCut", 0.0);
218 
220  trackQualityCuts.add<double>("chi2MaxPt", 10.)->setComment("max pT used to determine the pT-dependent chi2 cut");
221  trackQualityCuts.add<std::vector<double>>("chi2Coeff", {0.9, 1.8})->setComment("chi2 at 1GeV and at ptMax above");
222  trackQualityCuts.add<double>("chi2Scale", 8.)
223  ->setComment(
224  "Factor to multiply the pT-dependent chi2 cut (currently: 8 for the broken line fit, ?? for the Riemann "
225  "fit)");
226  trackQualityCuts.add<double>("tripletMinPt", 0.0)->setComment("Min pT for triplets, in GeV");
227  trackQualityCuts.add<double>("tripletMaxTip", 0.1)->setComment("Max |Tip| for triplets, in cm");
228  trackQualityCuts.add<double>("tripletMaxZip", 6.)->setComment("Max |Zip| for triplets, in cm");
229  trackQualityCuts.add<double>("quadrupletMinPt", 0.0)->setComment("Min pT for quadruplets, in GeV");
230  trackQualityCuts.add<double>("quadrupletMaxTip", 0.5)->setComment("Max |Tip| for quadruplets, in cm");
231  trackQualityCuts.add<double>("quadrupletMaxZip", 6.)->setComment("Max |Zip| for quadruplets, in cm");
232 
233  desc.add<std::vector<int>>(
234  "phiCuts", std::vector<int>(std::begin(phase1PixelTopology::phicuts), std::end(phase1PixelTopology::phicuts)))
235  ->setComment("Cuts in phi for cells");
236 
237  desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
238  ->setComment(
239  "Quality cuts based on the results of the track fit:\n - apply a pT-dependent chi2 cut;\n - apply \"region "
240  "cuts\" based on the fit results (pT, Tip, Zip).");
241 }
242 
243 template <>
245  fillDescriptionsCommon(desc);
246 
247  desc.add<bool>("idealConditions", false);
248  desc.add<bool>("includeFarForwards", true);
249  desc.add<bool>("includeJumpingForwardDoublets", true);
250  desc.add<double>("z0Cut", 7.5);
251  desc.add<double>("ptCut", 0.85);
252 
254  trackQualityCuts.add<double>("maxChi2", 5.)->setComment("Max normalized chi2");
255  trackQualityCuts.add<double>("minPt", 0.5)->setComment("Min pT in GeV");
256  trackQualityCuts.add<double>("maxTip", 0.3)->setComment("Max |Tip| in cm");
257  trackQualityCuts.add<double>("maxZip", 12.)->setComment("Max |Zip|, in cm");
258 
259  desc.add<std::vector<int>>(
260  "phiCuts", std::vector<int>(std::begin(phase2PixelTopology::phicuts), std::end(phase2PixelTopology::phicuts)))
261  ->setComment("Cuts in phi for cells");
262 
263  desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
264  ->setComment(
265  "Quality cuts based on the results of the track fit:\n - apply cuts based on the fit results (pT, Tip, "
266  "Zip).");
267 }
268 
269 template <typename TrackerTraits>
271  // 87 cm/GeV = 1/(3.8T * 0.3)
272  // take less than radius given by the hardPtCut and reject everything below
273  // auto hardCurvCut = 1.f/(0.35 * 87.f);
274  desc.add<double>("ptmin", 0.9)->setComment("Cut on minimum pt");
275  desc.add<double>("CAThetaCutBarrel", 0.002)->setComment("Cut on RZ alignement for Barrel");
276  desc.add<double>("CAThetaCutForward", 0.003)->setComment("Cut on RZ alignment for Forward");
277  desc.add<double>("hardCurvCut", 1. / (0.35 * 87.))->setComment("Cut on minimum curvature");
278  desc.add<double>("dcaCutInnerTriplet", 0.15)->setComment("Cut on origin radius when the inner hit is on BPix1");
279  desc.add<double>("dcaCutOuterTriplet", 0.25)->setComment("Cut on origin radius when the outer hit is on BPix1");
280  desc.add<bool>("earlyFishbone", true);
281  desc.add<bool>("lateFishbone", false);
282  desc.add<bool>("fillStatistics", false);
283  desc.add<unsigned int>("minHitsPerNtuplet", 4);
284  desc.add<unsigned int>("maxNumberOfDoublets", TrackerTraits::maxNumberOfDoublets);
285  desc.add<unsigned int>("minHitsForSharingCut", 10)
286  ->setComment("Maximum number of hits in a tuple to clean also if the shared hit is on bpx1");
287 
288  desc.add<bool>("fitNas4", false)->setComment("fit only 4 hits out of N");
289  desc.add<bool>("doClusterCut", true);
290  desc.add<bool>("doZ0Cut", true);
291  desc.add<bool>("doPtCut", true);
292  desc.add<bool>("useRiemannFit", false)->setComment("true for Riemann, false for BrokenLine");
293  desc.add<bool>("doSharedHitCut", true)->setComment("Sharing hit nTuples cleaning");
294  desc.add<bool>("dupPassThrough", false)->setComment("Do not reject duplicate");
295  desc.add<bool>("useSimpleTripletCleaner", true)->setComment("use alternate implementation");
296 }
297 
298 template <typename TrackerTraits>
300  if (m_params.onGPU_) {
301  // allocate pinned host memory only if CUDA is available
303  if (cuda and cuda->enabled()) {
304  cudaCheck(cudaMalloc(&m_counters, sizeof(Counters)));
305  cudaCheck(cudaMemset(m_counters, 0, sizeof(Counters)));
306  }
307  } else {
308  m_counters = new Counters();
309  memset(m_counters, 0, sizeof(Counters));
310  }
311 }
312 
313 template <typename TrackerTraits>
315  if (m_params.onGPU_) {
316  // print the gpu statistics and free pinned host memory only if CUDA is available
318  if (cuda and cuda->enabled()) {
319  if (m_params.doStats_) {
320  // crash on multi-gpu processes
322  }
323  cudaFree(m_counters);
324  }
325  } else {
326  if (m_params.doStats_) {
328  }
329  delete m_counters;
330  }
331 }
332 
333 template <typename TrackerTraits>
335  HitsOnDevice const& hits_d, float bfield, cudaStream_t stream) const {
339 
341 
342  GPUKernels kernels(m_params);
343  kernels.setCounters(m_counters);
344  kernels.allocateOnGPU(hits_d.nHits(), stream);
345 
346  kernels.buildDoublets(hits_d.view(), hits_d.offsetBPIX2(), stream);
347 
348  kernels.launchKernels(hits_d.view(), tracks.view(), stream);
349 
350  HelixFitOnGPU fitter(bfield, m_params.fitNas4_);
351  fitter.allocateOnGPU(kernels.tupleMultiplicity(), tracks.view());
352  if (m_params.useRiemannFit_) {
353  fitter.launchRiemannKernels(hits_d.view(), hits_d.nHits(), TrackerTraits::maxNumberOfQuadruplets, stream);
354  } else {
355  fitter.launchBrokenLineKernels(hits_d.view(), hits_d.nHits(), TrackerTraits::maxNumberOfQuadruplets, stream);
356  }
357  kernels.classifyTuples(hits_d.view(), tracks.view(), stream);
358 #ifdef GPU_DEBUG
359  cudaDeviceSynchronize();
360  cudaCheck(cudaGetLastError());
361  std::cout << "finished building pixel tracks on GPU" << std::endl;
362 #endif
363 
364  return tracks;
365 }
366 
367 template <typename TrackerTraits>
369  float bfield) const {
373 
375 
376  CPUKernels kernels(m_params);
377  kernels.setCounters(m_counters);
378  kernels.allocateOnGPU(hits_h.nHits(), nullptr);
379 
380  kernels.buildDoublets(hits_h.view(), hits_h.offsetBPIX2(), nullptr);
381  kernels.launchKernels(hits_h.view(), tracks.view(), nullptr);
382 
383  if (0 == hits_h.nHits())
384  return tracks;
385 
386  // now fit
387  HelixFitOnGPU fitter(bfield, m_params.fitNas4_);
388  fitter.allocateOnGPU(kernels.tupleMultiplicity(), tracks.view());
389 
390  if (m_params.useRiemannFit_) {
391  fitter.launchRiemannKernelsOnCPU(hits_h.view(), hits_h.nHits(), TrackerTraits::maxNumberOfQuadruplets);
392  } else {
393  fitter.launchBrokenLineKernelsOnCPU(hits_h.view(), hits_h.nHits(), TrackerTraits::maxNumberOfQuadruplets);
394  }
395 
396  kernels.classifyTuples(hits_h.view(), tracks.view(), nullptr);
397 
398 #ifdef GPU_DEBUG
399  std::cout << "finished building pixel tracks on CPU" << std::endl;
400 #endif
401 
402  // check that the fixed-size SoA does not overflow
403  auto maxTracks = tracks.view().metadata().size();
404  auto nTracks = tracks.view().nTracks();
406  if (nTracks == maxTracks - 1) {
407  edm::LogWarning("PixelTracks") << "Unsorted reconstructed pixel tracks truncated to " << maxTracks - 1
408  << " candidates";
409  }
410 
411  return tracks;
412 }
413 
static void fillDescriptionsCommon(edm::ParameterSetDescription &desc)
static void printCounters(Counters const *counters)
TrackSoAHost makeTuples(HitsOnHost const &hits_d, float bfield) const
CAHitNtupletGeneratorOnGPU(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
constexpr int16_t phicuts[nPairs]
typename std::enable_if< std::is_base_of< Phase2, T >::value >::type isPhase2Topology
caHitNtupletGenerator::Counters Counters
static void fillDescriptions(edm::ParameterSetDescription &desc)
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t stream
assert(be >=bs)
uint32_t offsetBPIX2() const
void launchRiemannKernelsOnCPU(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples)
TrackSoADevice makeTuplesAsync(HitsOnDevice const &hits_d, float bfield, cudaStream_t stream) const
void launchRiemannKernels(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples, cudaStream_t cudaStream)
typename std::enable_if< std::is_base_of< Phase1, T >::value >::type isPhase1Topology
constexpr int16_t phicuts[nPairs]
void launchBrokenLineKernelsOnCPU(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples)
void allocateOnGPU(TupleMultiplicity const *tupleMultiplicity, OutputSoAView &helix_fit_results)
Definition: HelixFitOnGPU.cc:5
ParameterSet const & getParameterSet(ParameterSetID const &id)
maxTracks
Definition: DMR_cfg.py:158
Square< F >::type sqr(const F &f)
Definition: Square.h:14
float x
#define cudaCheck(ARG,...)
Definition: cudaCheck.h:69
Log< level::Warning, false > LogWarning
long double T
void launchBrokenLineKernels(const HitConstView &hv, uint32_t nhits, uint32_t maxNumberOfTuples, cudaStream_t cudaStream)
static void printCounters(Counters const *counters)