CMS 3D CMS Logo

CAHitNtupletGenerator.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 
19 
20 #include "CAHitNtupletGenerator.h"
22 #include "CAPixelDoublets.h"
23 #include "CAPixelDoubletsAlgos.h"
24 
26  namespace {
27 
28  using namespace caHitNtupletGenerator;
29  using namespace caPixelDoublets;
30  using namespace pixelTopology;
31  using namespace pixelTrack;
32 
33  template <typename T>
34  T sqr(T x) {
35  return x * x;
36  }
37 
38  //Common Params
39  void fillDescriptionsCommon(edm::ParameterSetDescription& desc) {
40  // 87 cm/GeV = 1/(3.8T * 0.3)
41  // take less than radius given by the hardPtCut and reject everything below
42  // auto hardCurvCut = 1.f/(0.35 * 87.f);
43  desc.add<double>("ptmin", 0.9f)->setComment("Cut on minimum pt");
44  desc.add<double>("CAThetaCutBarrel", 0.002f)->setComment("Cut on RZ alignement for Barrel");
45  desc.add<double>("CAThetaCutForward", 0.003f)->setComment("Cut on RZ alignment for Forward");
46  desc.add<double>("hardCurvCut", 1.f / (0.35 * 87.f))
47  ->setComment("Cut on minimum curvature, used in DCA ntuplet selection");
48  desc.add<double>("dcaCutInnerTriplet", 0.15f)->setComment("Cut on origin radius when the inner hit is on BPix1");
49  desc.add<double>("dcaCutOuterTriplet", 0.25f)->setComment("Cut on origin radius when the outer hit is on BPix1");
50  desc.add<bool>("earlyFishbone", true);
51  desc.add<bool>("lateFishbone", false);
52  desc.add<bool>("fillStatistics", false);
53  desc.add<unsigned int>("minHitsPerNtuplet", 4);
54  desc.add<unsigned int>("minHitsForSharingCut", 10)
55  ->setComment("Maximum number of hits in a tuple to clean also if the shared hit is on bpx1");
56 
57  desc.add<bool>("fitNas4", false)->setComment("fit only 4 hits out of N");
58  desc.add<bool>("doClusterCut", true);
59  desc.add<bool>("doZ0Cut", true);
60  desc.add<bool>("doPtCut", true);
61  desc.add<bool>("useRiemannFit", false)->setComment("true for Riemann, false for BrokenLine");
62  desc.add<bool>("doSharedHitCut", true)->setComment("Sharing hit nTuples cleaning");
63  desc.add<bool>("dupPassThrough", false)->setComment("Do not reject duplicate");
64  desc.add<bool>("useSimpleTripletCleaner", true)->setComment("use alternate implementation");
65  }
66 
67  AlgoParams makeCommonParams(edm::ParameterSet const& cfg) {
68  return AlgoParams({cfg.getParameter<unsigned int>("minHitsForSharingCut"),
69  cfg.getParameter<bool>("useRiemannFit"),
70  cfg.getParameter<bool>("fitNas4"),
71  cfg.getParameter<bool>("includeJumpingForwardDoublets"),
72  cfg.getParameter<bool>("earlyFishbone"),
73  cfg.getParameter<bool>("lateFishbone"),
74  cfg.getParameter<bool>("fillStatistics"),
75  cfg.getParameter<bool>("doSharedHitCut"),
76  cfg.getParameter<bool>("dupPassThrough"),
77  cfg.getParameter<bool>("useSimpleTripletCleaner")});
78  }
79 
80  //This is needed to have the partial specialization for isPhase1Topology/isPhase2Topology
81  template <typename TrackerTraits, typename Enable = void>
82  struct TopologyCuts {};
83 
84  template <typename TrackerTraits>
85  struct TopologyCuts<TrackerTraits, isPhase1Topology<TrackerTraits>> {
86  static constexpr CAParamsT<TrackerTraits> makeCACuts(edm::ParameterSet const& cfg) {
87  return CAParamsT<TrackerTraits>{{cfg.getParameter<unsigned int>("maxNumberOfDoublets"),
88  cfg.getParameter<unsigned int>("minHitsPerNtuplet"),
89  (float)cfg.getParameter<double>("ptmin"),
90  (float)cfg.getParameter<double>("CAThetaCutBarrel"),
91  (float)cfg.getParameter<double>("CAThetaCutForward"),
92  (float)cfg.getParameter<double>("hardCurvCut"),
93  (float)cfg.getParameter<double>("dcaCutInnerTriplet"),
94  (float)cfg.getParameter<double>("dcaCutOuterTriplet")}};
95  };
96 
97  static constexpr ::pixelTrack::QualityCutsT<TrackerTraits> makeQualityCuts(edm::ParameterSet const& pset) {
98  auto coeff = pset.getParameter<std::array<double, 2>>("chi2Coeff");
99  auto ptMax = pset.getParameter<double>("chi2MaxPt");
100 
101  coeff[1] = (coeff[1] - coeff[0]) / log2(ptMax);
102  return ::pixelTrack::QualityCutsT<TrackerTraits>{// polynomial coefficients for the pT-dependent chi2 cut
103  {(float)coeff[0], (float)coeff[1], 0.f, 0.f},
104  // max pT used to determine the chi2 cut
105  (float)ptMax,
106  // chi2 scale factor: 8 for broken line fit, ?? for Riemann fit
107  (float)pset.getParameter<double>("chi2Scale"),
108  // regional cuts for triplets
109  {(float)pset.getParameter<double>("tripletMaxTip"),
110  (float)pset.getParameter<double>("tripletMinPt"),
111  (float)pset.getParameter<double>("tripletMaxZip")},
112  // regional cuts for quadruplets
113  {(float)pset.getParameter<double>("quadrupletMaxTip"),
114  (float)pset.getParameter<double>("quadrupletMinPt"),
115  (float)pset.getParameter<double>("quadrupletMaxZip")}};
116  }
117  };
118 
119  template <typename TrackerTraits>
120  struct TopologyCuts<TrackerTraits, isPhase2Topology<TrackerTraits>> {
121  static constexpr CAParamsT<TrackerTraits> makeCACuts(edm::ParameterSet const& cfg) {
122  return CAParamsT<TrackerTraits>{{cfg.getParameter<unsigned int>("maxNumberOfDoublets"),
123  cfg.getParameter<unsigned int>("minHitsPerNtuplet"),
124  (float)cfg.getParameter<double>("ptmin"),
125  (float)cfg.getParameter<double>("CAThetaCutBarrel"),
126  (float)cfg.getParameter<double>("CAThetaCutForward"),
127  (float)cfg.getParameter<double>("hardCurvCut"),
128  (float)cfg.getParameter<double>("dcaCutInnerTriplet"),
129  (float)cfg.getParameter<double>("dcaCutOuterTriplet")},
130  {(bool)cfg.getParameter<bool>("includeFarForwards")}};
131  }
132 
133  static constexpr ::pixelTrack::QualityCutsT<TrackerTraits> makeQualityCuts(edm::ParameterSet const& pset) {
134  return ::pixelTrack::QualityCutsT<TrackerTraits>{
135  static_cast<float>(pset.getParameter<double>("maxChi2")),
136  static_cast<float>(pset.getParameter<double>("minPt")),
137  static_cast<float>(pset.getParameter<double>("maxTip")),
138  static_cast<float>(pset.getParameter<double>("maxZip")),
139  };
140  }
141  };
142 
143  //Cell Cuts, as they are the cuts have the same logic for Phase2 and Phase1
144  //keeping them separate would allow further differentiation in the future
145  //moving them to TopologyCuts and using the same syntax
146  template <typename TrackerTraits>
147  CellCutsT<TrackerTraits> makeCellCuts(edm::ParameterSet const& cfg) {
148  return CellCutsT<TrackerTraits>{cfg.getParameter<bool>("doClusterCut"),
149  cfg.getParameter<bool>("doZ0Cut"),
150  cfg.getParameter<bool>("doPtCut"),
151  cfg.getParameter<bool>("idealConditions"),
152  (float)cfg.getParameter<double>("cellZ0Cut"),
153  (float)cfg.getParameter<double>("cellPtCut"),
154  cfg.getParameter<std::vector<int>>("phiCuts")};
155  }
156 
157  } // namespace
158 
159  using namespace std;
160 
161  template <typename TrackerTraits>
163  : m_params(makeCommonParams(cfg),
164  makeCellCuts<TrackerTraits>(cfg),
165  TopologyCuts<TrackerTraits>::makeQualityCuts(cfg.getParameterSet("trackQualityCuts")),
166  TopologyCuts<TrackerTraits>::makeCACuts(cfg)) {
167 #ifdef DUMP_GPU_TK_TUPLES
168  printf("TK: %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s\n",
169  "tid",
170  "qual",
171  "nh",
172  "nl",
173  "charge",
174  "pt",
175  "eta",
176  "phi",
177  "tip",
178  "zip",
179  "chi2",
180  "h1",
181  "h2",
182  "h3",
183  "h4",
184  "h5",
185  "hn");
186 #endif
187  }
188 
189  template <typename TrackerTraits>
191  static_assert(sizeof(TrackerTraits) == 0,
192  "Note: this fillPSetDescription is a dummy one. Please specialise it for the correct version of "
193  "CAHitNtupletGenerator<TrackerTraits>.");
194  }
195 
196  template <>
198  fillDescriptionsCommon(desc);
199 
200  desc.add<unsigned int>("maxNumberOfDoublets", pixelTopology::Phase1::maxNumberOfDoublets);
201  desc.add<bool>("idealConditions", true);
202  desc.add<bool>("includeJumpingForwardDoublets", false);
203  desc.add<double>("cellZ0Cut", 12.0);
204  desc.add<double>("cellPtCut", 0.5);
205 
207  trackQualityCuts.add<double>("chi2MaxPt", 10.)->setComment("max pT used to determine the pT-dependent chi2 cut");
208  trackQualityCuts.add<std::vector<double>>("chi2Coeff", {0.9, 1.8})->setComment("chi2 at 1GeV and at ptMax above");
209  trackQualityCuts.add<double>("chi2Scale", 8.)
210  ->setComment(
211  "Factor to multiply the pT-dependent chi2 cut (currently: 8 for the broken line fit, ?? for the Riemann "
212  "fit)");
213  trackQualityCuts.add<double>("tripletMinPt", 0.5)->setComment("Min pT for triplets, in GeV");
214  trackQualityCuts.add<double>("tripletMaxTip", 0.3)->setComment("Max |Tip| for triplets, in cm");
215  trackQualityCuts.add<double>("tripletMaxZip", 12.)->setComment("Max |Zip| for triplets, in cm");
216  trackQualityCuts.add<double>("quadrupletMinPt", 0.3)->setComment("Min pT for quadruplets, in GeV");
217  trackQualityCuts.add<double>("quadrupletMaxTip", 0.5)->setComment("Max |Tip| for quadruplets, in cm");
218  trackQualityCuts.add<double>("quadrupletMaxZip", 12.)->setComment("Max |Zip| for quadruplets, in cm");
219  desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
220  ->setComment(
221  "Quality cuts based on the results of the track fit:\n - apply a pT-dependent chi2 cut;\n - apply "
222  "\"region "
223  "cuts\" based on the fit results (pT, Tip, Zip).");
224 
225  desc.add<std::vector<int>>(
226  "phiCuts",
227  std::vector<int>(std::begin(phase1PixelTopology::phicuts), std::end(phase1PixelTopology::phicuts)))
228  ->setComment("Cuts in phi for cells");
229  }
230 
231  template <>
233  fillDescriptionsCommon(desc);
234 
235  desc.add<unsigned int>("maxNumberOfDoublets", pixelTopology::HIonPhase1::maxNumberOfDoublets);
236  desc.add<bool>("idealConditions", false);
237  desc.add<bool>("includeJumpingForwardDoublets", false);
238  desc.add<double>("cellZ0Cut", 10.0);
239  desc.add<double>("cellPtCut", 0.0);
240 
242  trackQualityCuts.add<double>("chi2MaxPt", 10.)->setComment("max pT used to determine the pT-dependent chi2 cut");
243  trackQualityCuts.add<std::vector<double>>("chi2Coeff", {0.9, 1.8})->setComment("chi2 at 1GeV and at ptMax above");
244  trackQualityCuts.add<double>("chi2Scale", 8.)
245  ->setComment(
246  "Factor to multiply the pT-dependent chi2 cut (currently: 8 for the broken line fit, ?? for the Riemann "
247  "fit)");
248  trackQualityCuts.add<double>("tripletMinPt", 0.0)->setComment("Min pT for triplets, in GeV");
249  trackQualityCuts.add<double>("tripletMaxTip", 0.1)->setComment("Max |Tip| for triplets, in cm");
250  trackQualityCuts.add<double>("tripletMaxZip", 6.)->setComment("Max |Zip| for triplets, in cm");
251  trackQualityCuts.add<double>("quadrupletMinPt", 0.0)->setComment("Min pT for quadruplets, in GeV");
252  trackQualityCuts.add<double>("quadrupletMaxTip", 0.5)->setComment("Max |Tip| for quadruplets, in cm");
253  trackQualityCuts.add<double>("quadrupletMaxZip", 6.)->setComment("Max |Zip| for quadruplets, in cm");
254 
255  desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
256  ->setComment(
257  "Quality cuts based on the results of the track fit:\n - apply a pT-dependent chi2 cut;\n - apply "
258  "\"region "
259  "cuts\" based on the fit results (pT, Tip, Zip).");
260 
261  desc.add<std::vector<int>>(
262  "phiCuts",
263  std::vector<int>(std::begin(phase1PixelTopology::phicuts), std::end(phase1PixelTopology::phicuts)))
264  ->setComment("Cuts in phi for cells");
265  }
266 
267  template <>
269  fillDescriptionsCommon(desc);
270 
271  desc.add<unsigned int>("maxNumberOfDoublets", pixelTopology::Phase2::maxNumberOfDoublets);
272  desc.add<bool>("idealConditions", false);
273  desc.add<bool>("includeFarForwards", true);
274  desc.add<bool>("includeJumpingForwardDoublets", true);
275  desc.add<double>("cellZ0Cut", 7.5);
276  desc.add<double>("cellPtCut", 0.85);
277 
279  trackQualityCuts.add<double>("maxChi2", 5.)->setComment("Max normalized chi2");
280  trackQualityCuts.add<double>("minPt", 0.5)->setComment("Min pT in GeV");
281  trackQualityCuts.add<double>("maxTip", 0.3)->setComment("Max |Tip| in cm");
282  trackQualityCuts.add<double>("maxZip", 12.)->setComment("Max |Zip|, in cm");
283  desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
284  ->setComment(
285  "Quality cuts based on the results of the track fit:\n - apply cuts based on the fit results (pT, Tip, "
286  "Zip).");
287 
288  desc.add<std::vector<int>>(
289  "phiCuts",
290  std::vector<int>(std::begin(phase2PixelTopology::phicuts), std::end(phase2PixelTopology::phicuts)))
291  ->setComment("Cuts in phi for cells");
292  }
293 
294  template <typename TrackerTraits>
296  HitsOnDevice const& hits_d, ParamsOnDevice const* cpeParams, float bfield, Queue& queue) const {
300 
302 
303  GPUKernels kernels(m_params, hits_d.view().metadata().size(), queue);
304 
305  kernels.buildDoublets(hits_d.view(), queue);
306  kernels.launchKernels(hits_d.view(), tracks.view(), queue);
307 
308  HelixFit fitter(bfield, m_params.fitNas4_);
309  fitter.allocate(kernels.tupleMultiplicity(), tracks.view());
310  if (m_params.useRiemannFit_) {
311  fitter.launchRiemannKernels(
312  hits_d.view(), cpeParams, hits_d.view().metadata().size(), TrackerTraits::maxNumberOfQuadruplets, queue);
313  } else {
315  hits_d.view(), cpeParams, hits_d.view().metadata().size(), TrackerTraits::maxNumberOfQuadruplets, queue);
316  }
317  kernels.classifyTuples(hits_d.view(), tracks.view(), queue);
318 #ifdef GPU_DEBUG
320  std::cout << "finished building pixel tracks on GPU" << std::endl;
321 #endif
322 
323  return tracks;
324  }
325 
329 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
constexpr int16_t phicuts[nPairs]
typename std::enable_if< std::is_base_of< Phase2, T >::value >::type isPhase2Topology
void launchBrokenLineKernels(const HitConstView &hv, ParamsOnDevice const *cpeParams, uint32_t nhits, uint32_t maxNumberOfTuples, Queue &queue)
TrackingRecHitsSoACollection< TrackerTraits > HitsOnDevice
typename std::enable_if< std::is_base_of< Phase1, T >::value >::type isPhase1Topology
TkSoADevice makeTuplesAsync(HitsOnDevice const &hits_d, ParamsOnDevice const *cpeParams, float bfield, Queue &queue) const
constexpr int16_t phicuts[nPairs]
static constexpr uint32_t maxNumberOfDoublets
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, TracksHost< TrackerTraits >, TracksDevice< TrackerTraits, Device > > TracksSoACollection
ParameterSet const & getParameterSet(ParameterSetID const &id)
static constexpr uint32_t maxNumberOfDoublets
void launchRiemannKernels(const HitConstView &hv, ParamsOnDevice const *cpeParams, uint32_t nhits, uint32_t maxNumberOfTuples, Queue &queue)
static constexpr uint32_t maxNumberOfDoublets
Square< F >::type sqr(const F &f)
Definition: Square.h:14
static void fillPSetDescription(edm::ParameterSetDescription &desc)
float x
void allocate(TupleMultiplicity const *tupleMultiplicity, OutputSoAView &helix_fit_results)
Definition: HelixFit.cc:6
long double T