File indexing completed on 2025-07-03 00:42:41
0001
0002
0003
0004 #include <array>
0005 #include <cassert>
0006 #include <functional>
0007 #include <vector>
0008
0009 #include <alpaka/alpaka.hpp>
0010
0011 #include "DataFormats/TrackSoA/interface/TracksDevice.h"
0012 #include "DataFormats/TrackSoA/interface/TracksHost.h"
0013 #include "DataFormats/TrackSoA/interface/alpaka/TracksSoACollection.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0016 #include "FWCore/Utilities/interface/Exception.h"
0017 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0018
0019 #include "CAHitNtupletGenerator.h"
0020 #include "CAHitNtupletGeneratorKernels.h"
0021 #include "CAPixelDoublets.h"
0022 #include "CAPixelDoubletsAlgos.h"
0023
0024 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0025 namespace {
0026
0027 using namespace caHitNtupletGenerator;
0028 using namespace caPixelDoublets;
0029 using namespace caStructures;
0030 using namespace pixelTopology;
0031 using namespace pixelTrack;
0032
0033 template <typename T>
0034 T sqr(T x) {
0035 return x * x;
0036 }
0037
0038
0039 void fillDescriptionsCommon(edm::ParameterSetDescription& desc) {
0040 desc.add<double>("cellZ0Cut", 12.0f)->setComment("Z0 cut for cells");
0041 desc.add<double>("cellPtCut", 0.5f)->setComment("Preliminary pT cut at cell building level.");
0042
0043
0044 desc.add<double>("dzdrFact", 8.0f * 0.0285f / 0.015f);
0045 desc.add<int>("minYsizeB1", 1)
0046 ->setComment("Cut on inner hit cluster size (in Y) for barrel-forward cells. Barrel 1 cut.");
0047 desc.add<int>("minYsizeB2", 1)
0048 ->setComment("Cut on inner hit cluster size (in Y) for barrel-forward cells. Barrel 2 cut.");
0049 desc.add<int>("maxDYsize12", 28)
0050 ->setComment("Cut on cluster size differences (in Y) for barrel-forward cells. Barrel 1-2 cells.");
0051 desc.add<int>("maxDYsize", 20)
0052 ->setComment("Cut on cluster size differences (in Y) for barrel-forward cells. Other barrel cells.");
0053 desc.add<int>("maxDYPred", 20)
0054 ->setComment("Cut on cluster size differences (in Y) for barrel-forward cells. Barrel-forward cells.");
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068
0069
0070 desc.add<std::string>("maxNumberOfDoublets", std::to_string(pixelTopology::Phase1::maxNumberOfDoublets))
0071 ->setComment(
0072 "Max nummber of doublets (cells) as a string. The string will be parsed to a TFormula, depending on "
0073 "nHits (labeled 'x'), \n and evaluated for each event. May also be a constant.");
0074 desc.add<std::string>("maxNumberOfTuples", std::to_string(pixelTopology::Phase1::maxNumberOfTuples))
0075 ->setComment("Max nummber of tuples as a string. Same behavior as maxNumberOfDoublets.");
0076 desc.add<double>("avgHitsPerTrack", 5.0f)->setComment("Number of hits per track. Average per track.");
0077 desc.add<double>("avgCellsPerHit", 25.0f)
0078 ->setComment("Number of cells for which an hit is the outer hit. Average per hit.");
0079 desc.add<double>("avgCellsPerCell", 2.0f)
0080 ->setComment("Number of cells connected to another cell. Average per cell.");
0081 desc.add<double>("avgTracksPerCell", 1.0f)
0082 ->setComment("Number of tracks to which a cell belongs. Average per cell.");
0083
0084
0085 desc.add<double>("ptmin", 0.9f)->setComment("Cut on minimum pt");
0086
0087
0088
0089 desc.add<double>("hardCurvCut", 1.f / (0.35 * 87.f))
0090 ->setComment("Cut on minimum curvature, used in DCA ntuplet selection");
0091
0092 desc.add<bool>("earlyFishbone", true);
0093 desc.add<bool>("lateFishbone", false);
0094 desc.add<bool>("fillStatistics", false);
0095 desc.add<unsigned int>("minHitsPerNtuplet", 4);
0096 desc.add<unsigned int>("minHitsForSharingCut", 10)
0097 ->setComment("Maximum number of hits in a tuple to clean also if the shared hit is on bpx1");
0098
0099 desc.add<bool>("fitNas4", false)->setComment("fit only 4 hits out of N");
0100
0101 desc.add<bool>("useRiemannFit", false)->setComment("true for Riemann, false for BrokenLine");
0102 desc.add<bool>("doSharedHitCut", true)->setComment("Sharing hit nTuples cleaning");
0103 desc.add<bool>("dupPassThrough", false)->setComment("Do not reject duplicate");
0104 desc.add<bool>("useSimpleTripletCleaner", true)->setComment("use alternate implementation");
0105 }
0106
0107 AlgoParams makeCommonParams(edm::ParameterSet const& cfg) {
0108 return AlgoParams({
0109
0110
0111 (float)cfg.getParameter<double>("avgHitsPerTrack"),
0112 (float)cfg.getParameter<double>("avgCellsPerHit"),
0113 (float)cfg.getParameter<double>("avgCellsPerCell"),
0114 (float)cfg.getParameter<double>("avgTracksPerCell"),
0115
0116
0117 (uint16_t)cfg.getParameter<unsigned int>("minHitsPerNtuplet"),
0118 (uint16_t)cfg.getParameter<unsigned int>("minHitsForSharingCut"),
0119 (float)cfg.getParameter<double>("ptmin"),
0120 (float)cfg.getParameter<double>("hardCurvCut"),
0121 (float)cfg.getParameter<double>("cellZ0Cut"),
0122 (float)cfg.getParameter<double>("cellPtCut"),
0123
0124
0125 (float)cfg.getParameter<double>("dzdrFact"),
0126 (int16_t)cfg.getParameter<int>("minYsizeB1"),
0127 (int16_t)cfg.getParameter<int>("minYsizeB2"),
0128 (int16_t)cfg.getParameter<int>("maxDYsize12"),
0129 (int16_t)cfg.getParameter<int>("maxDYsize"),
0130 (int16_t)cfg.getParameter<int>("maxDYPred"),
0131
0132
0133 cfg.getParameter<bool>("useRiemannFit"),
0134 cfg.getParameter<bool>("fitNas4"),
0135 cfg.getParameter<bool>("earlyFishbone"),
0136 cfg.getParameter<bool>("lateFishbone"),
0137 cfg.getParameter<bool>("fillStatistics"),
0138 cfg.getParameter<bool>("doSharedHitCut"),
0139 cfg.getParameter<bool>("dupPassThrough"),
0140 cfg.getParameter<bool>("useSimpleTripletCleaner")});
0141 }
0142
0143
0144 template <typename TrackerTraits, typename Enable = void>
0145 struct TopologyCuts {};
0146
0147 template <typename TrackerTraits>
0148 struct TopologyCuts<TrackerTraits, isPhase1Topology<TrackerTraits>> {
0149 static constexpr ::pixelTrack::QualityCutsT<TrackerTraits> makeQualityCuts(edm::ParameterSet const& pset) {
0150 auto coeff = pset.getParameter<std::array<double, 2>>("chi2Coeff");
0151 auto ptMax = pset.getParameter<double>("chi2MaxPt");
0152
0153 coeff[1] = (coeff[1] - coeff[0]) / log2(ptMax);
0154 return ::pixelTrack::QualityCutsT<TrackerTraits>{
0155 {(float)coeff[0], (float)coeff[1], 0.f, 0.f},
0156
0157 (float)ptMax,
0158
0159 (float)pset.getParameter<double>("chi2Scale"),
0160
0161 {(float)pset.getParameter<double>("tripletMaxTip"),
0162 (float)pset.getParameter<double>("tripletMinPt"),
0163 (float)pset.getParameter<double>("tripletMaxZip")},
0164
0165 {(float)pset.getParameter<double>("quadrupletMaxTip"),
0166 (float)pset.getParameter<double>("quadrupletMinPt"),
0167 (float)pset.getParameter<double>("quadrupletMaxZip")}};
0168 }
0169 };
0170
0171 template <typename TrackerTraits>
0172 struct TopologyCuts<TrackerTraits, isPhase2Topology<TrackerTraits>> {
0173 static constexpr ::pixelTrack::QualityCutsT<TrackerTraits> makeQualityCuts(edm::ParameterSet const& pset) {
0174 return ::pixelTrack::QualityCutsT<TrackerTraits>{
0175 static_cast<float>(pset.getParameter<double>("maxChi2")),
0176 static_cast<float>(pset.getParameter<double>("minPt")),
0177 static_cast<float>(pset.getParameter<double>("maxTip")),
0178 static_cast<float>(pset.getParameter<double>("maxZip")),
0179 };
0180 }
0181 };
0182
0183 }
0184
0185 using namespace std;
0186
0187 template <typename TrackerTraits>
0188 CAHitNtupletGenerator<TrackerTraits>::CAHitNtupletGenerator(const edm::ParameterSet& cfg)
0189 : m_params(makeCommonParams(cfg),
0190 TopologyCuts<TrackerTraits>::makeQualityCuts(cfg.getParameterSet("trackQualityCuts"))) {
0191 #ifdef DUMP_GPU_TK_TUPLES
0192 printf("TK: %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s %s\n",
0193 "tid",
0194 "qual",
0195 "nh",
0196 "nl",
0197 "charge",
0198 "pt",
0199 "eta",
0200 "phi",
0201 "tip",
0202 "zip",
0203 "chi2",
0204 "h1",
0205 "h2",
0206 "h3",
0207 "h4",
0208 "h5",
0209 "hn");
0210 #endif
0211 }
0212
0213 template <typename TrackerTraits>
0214 void CAHitNtupletGenerator<TrackerTraits>::fillPSetDescription(edm::ParameterSetDescription& desc) {
0215 static_assert(sizeof(TrackerTraits) == 0,
0216 "Note: this fillPSetDescription is a dummy one. Please specialise it for the correct version of "
0217 "CAHitNtupletGenerator<TrackerTraits>.");
0218 }
0219
0220 template <>
0221 void CAHitNtupletGenerator<pixelTopology::Phase1>::fillPSetDescription(edm::ParameterSetDescription& desc) {
0222 fillDescriptionsCommon(desc);
0223
0224 edm::ParameterSetDescription trackQualityCuts;
0225 trackQualityCuts.add<double>("chi2MaxPt", 10.)->setComment("max pT used to determine the pT-dependent chi2 cut");
0226 trackQualityCuts.add<std::vector<double>>("chi2Coeff", {0.9, 1.8})->setComment("chi2 at 1GeV and at ptMax above");
0227 trackQualityCuts.add<double>("chi2Scale", 8.)
0228 ->setComment(
0229 "Factor to multiply the pT-dependent chi2 cut (currently: 8 for the broken line fit, ?? for the Riemann "
0230 "fit)");
0231 trackQualityCuts.add<double>("tripletMinPt", 0.5)->setComment("Min pT for triplets, in GeV");
0232 trackQualityCuts.add<double>("tripletMaxTip", 0.3)->setComment("Max |Tip| for triplets, in cm");
0233 trackQualityCuts.add<double>("tripletMaxZip", 12.)->setComment("Max |Zip| for triplets, in cm");
0234 trackQualityCuts.add<double>("quadrupletMinPt", 0.3)->setComment("Min pT for quadruplets, in GeV");
0235 trackQualityCuts.add<double>("quadrupletMaxTip", 0.5)->setComment("Max |Tip| for quadruplets, in cm");
0236 trackQualityCuts.add<double>("quadrupletMaxZip", 12.)->setComment("Max |Zip| for quadruplets, in cm");
0237 desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
0238 ->setComment(
0239 "Quality cuts based on the results of the track fit:\n - apply a pT-dependent chi2 cut;\n - apply "
0240 "\"region "
0241 "cuts\" based on the fit results (pT, Tip, Zip).");
0242
0243 edm::ParameterSetDescription geometryParams;
0244 using namespace phase1PixelTopology;
0245
0246 geometryParams
0247 .add<std::vector<double>>("caDCACuts",
0248 std::vector<double>(std::begin(dcaCuts), std::begin(dcaCuts) + numberOfLayers))
0249 ->setComment("Cut on RZ alignement. One per layer, the layer being the middle one for a triplet.");
0250 geometryParams
0251 .add<std::vector<double>>("caThetaCuts",
0252 std::vector<double>(std::begin(thetaCuts), std::begin(thetaCuts) + numberOfLayers))
0253 ->setComment("Cut on origin radius. One per layer, the layer being the innermost one for a triplet.");
0254 geometryParams.add<std::vector<unsigned int>>("startingPairs", {0u, 1u, 2u})
0255 ->setComment(
0256 "The list of the ids of pairs from which the CA ntuplets building may start.");
0257
0258 geometryParams
0259 .add<std::vector<unsigned int>>(
0260 "pairGraph",
0261 std::vector<unsigned int>(std::begin(layerPairs),
0262 std::begin(layerPairs) + (pixelTopology::Phase1::nPairsForQuadruplets * 2)))
0263 ->setComment("CA graph");
0264 geometryParams
0265 .add<std::vector<int>>(
0266 "phiCuts",
0267 std::vector<int>(std::begin(phicuts), std::begin(phicuts) + pixelTopology::Phase1::nPairsForQuadruplets))
0268 ->setComment("Cuts in phi for cells");
0269 geometryParams
0270 .add<std::vector<double>>(
0271 "minZ",
0272 std::vector<double>(std::begin(minz), std::begin(minz) + pixelTopology::Phase1::nPairsForQuadruplets))
0273 ->setComment("Cuts in min z (on inner hit) for cells");
0274 geometryParams
0275 .add<std::vector<double>>(
0276 "maxZ",
0277 std::vector<double>(std::begin(maxz), std::begin(maxz) + pixelTopology::Phase1::nPairsForQuadruplets))
0278 ->setComment("Cuts in max z (on inner hit) for cells");
0279 geometryParams
0280 .add<std::vector<double>>(
0281 "maxR",
0282 std::vector<double>(std::begin(maxr), std::begin(maxr) + pixelTopology::Phase1::nPairsForQuadruplets))
0283 ->setComment("Cuts in max r for cells");
0284
0285 desc.add<edm::ParameterSetDescription>("geometry", geometryParams)
0286 ->setComment(
0287 "Quality cuts based on the results of the track fit:\n - apply cuts based on the fit results (pT, Tip, "
0288 "Zip).");
0289 }
0290
0291 template <>
0292 void CAHitNtupletGenerator<pixelTopology::HIonPhase1>::fillPSetDescription(edm::ParameterSetDescription& desc) {
0293 fillDescriptionsCommon(desc);
0294
0295 edm::ParameterSetDescription trackQualityCuts;
0296 trackQualityCuts.add<double>("chi2MaxPt", 10.)->setComment("max pT used to determine the pT-dependent chi2 cut");
0297 trackQualityCuts.add<std::vector<double>>("chi2Coeff", {0.9, 1.8})->setComment("chi2 at 1GeV and at ptMax above");
0298 trackQualityCuts.add<double>("chi2Scale", 8.)
0299 ->setComment(
0300 "Factor to multiply the pT-dependent chi2 cut (currently: 8 for the broken line fit, ?? for the Riemann "
0301 "fit)");
0302 trackQualityCuts.add<double>("tripletMinPt", 0.0)->setComment("Min pT for triplets, in GeV");
0303 trackQualityCuts.add<double>("tripletMaxTip", 0.1)->setComment("Max |Tip| for triplets, in cm");
0304 trackQualityCuts.add<double>("tripletMaxZip", 6.)->setComment("Max |Zip| for triplets, in cm");
0305 trackQualityCuts.add<double>("quadrupletMinPt", 0.0)->setComment("Min pT for quadruplets, in GeV");
0306 trackQualityCuts.add<double>("quadrupletMaxTip", 0.5)->setComment("Max |Tip| for quadruplets, in cm");
0307 trackQualityCuts.add<double>("quadrupletMaxZip", 6.)->setComment("Max |Zip| for quadruplets, in cm");
0308
0309 desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
0310 ->setComment(
0311 "Quality cuts based on the results of the track fit:\n - apply a pT-dependent chi2 cut;\n - apply "
0312 "\"region "
0313 "cuts\" based on the fit results (pT, Tip, Zip).");
0314
0315 edm::ParameterSetDescription geometryParams;
0316 using namespace phase1PixelTopology;
0317
0318 geometryParams
0319 .add<std::vector<double>>("caDCACuts",
0320 std::vector<double>(std::begin(phase1HIonPixelTopology::dcaCuts),
0321 std::begin(phase1HIonPixelTopology::dcaCuts) + numberOfLayers))
0322 ->setComment("Cut on RZ alignement. One per layer, the layer being the middle one for a triplet.");
0323 geometryParams
0324 .add<std::vector<double>>("caThetaCuts",
0325 std::vector<double>(std::begin(phase1HIonPixelTopology::thetaCuts),
0326 std::begin(phase1HIonPixelTopology::thetaCuts) + numberOfLayers))
0327 ->setComment("Cut on origin radius. One per layer, the layer being the innermost one for a triplet.");
0328 geometryParams.add<std::vector<unsigned int>>("startingPairs", {0u, 1u, 2u})
0329 ->setComment(
0330 "The list of the ids of pairs from which the CA ntuplets building may start.");
0331
0332 geometryParams
0333 .add<std::vector<unsigned int>>(
0334 "pairGraph",
0335 std::vector<unsigned int>(std::begin(layerPairs),
0336 std::begin(layerPairs) + (pixelTopology::Phase1::nPairsForQuadruplets * 2)))
0337 ->setComment("CA graph");
0338 geometryParams
0339 .add<std::vector<int>>("phiCuts",
0340 std::vector<int>(std::begin(phase1HIonPixelTopology::phicuts),
0341 std::begin(phase1HIonPixelTopology::phicuts) +
0342 pixelTopology::Phase1::nPairsForQuadruplets))
0343 ->setComment("Cuts in phi for cells");
0344 geometryParams
0345 .add<std::vector<double>>(
0346 "minZ",
0347 std::vector<double>(std::begin(minz), std::begin(minz) + pixelTopology::Phase1::nPairsForQuadruplets))
0348 ->setComment("Cuts in min z (on inner hit) for cells");
0349 geometryParams
0350 .add<std::vector<double>>(
0351 "maxZ",
0352 std::vector<double>(std::begin(maxz), std::begin(maxz) + pixelTopology::Phase1::nPairsForQuadruplets))
0353 ->setComment("Cuts in max z (on inner hit) for cells");
0354 geometryParams
0355 .add<std::vector<double>>(
0356 "maxR",
0357 std::vector<double>(std::begin(maxr), std::begin(maxr) + pixelTopology::Phase1::nPairsForQuadruplets))
0358 ->setComment("Cuts in max r for cells");
0359
0360 desc.add<edm::ParameterSetDescription>("geometry", geometryParams)
0361 ->setComment(
0362 "Quality cuts based on the results of the track fit:\n - apply cuts based on the fit results (pT, Tip, "
0363 "Zip).");
0364 }
0365
0366 template <>
0367 void CAHitNtupletGenerator<pixelTopology::Phase2>::fillPSetDescription(edm::ParameterSetDescription& desc) {
0368 fillDescriptionsCommon(desc);
0369
0370 edm::ParameterSetDescription trackQualityCuts;
0371 trackQualityCuts.add<double>("maxChi2", 5.)->setComment("Max normalized chi2");
0372 trackQualityCuts.add<double>("minPt", 0.5)->setComment("Min pT in GeV");
0373 trackQualityCuts.add<double>("maxTip", 0.3)->setComment("Max |Tip| in cm");
0374 trackQualityCuts.add<double>("maxZip", 12.)->setComment("Max |Zip|, in cm");
0375 desc.add<edm::ParameterSetDescription>("trackQualityCuts", trackQualityCuts)
0376 ->setComment(
0377 "Quality cuts based on the results of the track fit:\n - apply cuts based on the fit results (pT, Tip, "
0378 "Zip).");
0379
0380 edm::ParameterSetDescription geometryParams;
0381 using namespace phase2PixelTopology;
0382
0383 geometryParams
0384 .add<std::vector<double>>("caDCACuts",
0385 std::vector<double>(std::begin(dcaCuts), std::begin(dcaCuts) + numberOfLayers))
0386 ->setComment("Cut on RZ alignement. One per layer, the layer being the middle one for a triplet.");
0387 geometryParams
0388 .add<std::vector<double>>("caThetaCuts",
0389 std::vector<double>(std::begin(thetaCuts), std::begin(thetaCuts) + numberOfLayers))
0390 ->setComment("Cut on origin radius. One per layer, the layer being the innermost one for a triplet.");
0391 geometryParams
0392 .add<std::vector<unsigned int>>("startingPairs",
0393 {0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16,
0394 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32})
0395 ->setComment(
0396 "The list of the ids of pairs from which the CA ntuplets building may start.");
0397
0398 geometryParams
0399 .add<std::vector<unsigned int>>(
0400 "pairGraph", std::vector<unsigned int>(std::begin(layerPairs), std::begin(layerPairs) + (nPairs * 2)))
0401 ->setComment("CA graph");
0402 geometryParams
0403 .add<std::vector<int>>("phiCuts", std::vector<int>(std::begin(phicuts), std::begin(phicuts) + nPairs))
0404 ->setComment("Cuts in phi for cells");
0405 geometryParams.add<std::vector<double>>("minZ", std::vector<double>(std::begin(minz), std::begin(minz) + nPairs))
0406 ->setComment("Cuts in min z (on inner hit) for cells");
0407 geometryParams.add<std::vector<double>>("maxZ", std::vector<double>(std::begin(maxz), std::begin(maxz) + nPairs))
0408 ->setComment("Cuts in max z (on inner hit) for cells");
0409 geometryParams.add<std::vector<double>>("maxR", std::vector<double>(std::begin(maxr), std::begin(maxr) + nPairs))
0410 ->setComment("Cuts in max r for cells");
0411
0412 desc.add<edm::ParameterSetDescription>("geometry", geometryParams)
0413 ->setComment(
0414 "Quality cuts based on the results of the track fit:\n - apply cuts based on the fit results (pT, Tip, "
0415 "Zip).");
0416 }
0417
0418 template <typename TrackerTraits>
0419 reco::TracksSoACollection CAHitNtupletGenerator<TrackerTraits>::makeTuplesAsync(HitsOnDevice const& hits_d,
0420 CAGeometryOnDevice const& geometry_d,
0421 float bfield,
0422 uint32_t nDoublets,
0423 uint32_t nTracks,
0424 Queue& queue) const {
0425 using HelixFit = HelixFit<TrackerTraits>;
0426 using GPUKernels = CAHitNtupletGeneratorKernels<TrackerTraits>;
0427 using TrackHitSoA = ::reco::TrackHitSoA;
0428 using HitContainer = caStructures::HitContainerT<TrackerTraits>;
0429
0430 const int32_t H = m_params.algoParams_.avgHitsPerTrack_;
0431
0432 reco::TracksSoACollection tracks({{int(nTracks), int(nTracks * H)}}, queue);
0433
0434
0435 if (hits_d.view().metadata().size() < 2) {
0436 const auto device = alpaka::getDev(queue);
0437 auto ntracks_d = cms::alpakatools::make_device_view(device, tracks.view().nTracks());
0438 alpaka::memset(queue, ntracks_d, 0);
0439 return tracks;
0440 }
0441 GPUKernels kernels(
0442 m_params, hits_d.nHits(), hits_d.offsetBPIX2(), nDoublets, nTracks, geometry_d.view().metadata().size(), queue);
0443
0444 kernels.prepareHits(hits_d.view(), hits_d.view<::reco::HitModuleSoA>(), geometry_d.view(), queue);
0445 kernels.buildDoublets(hits_d.view(),
0446 geometry_d.view<::reco::CAGraphSoA>(),
0447 geometry_d.view<::reco::CALayersSoA>(),
0448 hits_d.offsetBPIX2(),
0449 queue);
0450 kernels.launchKernels(hits_d.view(),
0451 hits_d.offsetBPIX2(),
0452 geometry_d.view().metadata().size(),
0453 tracks.view(),
0454 tracks.view<TrackHitSoA>(),
0455 geometry_d.view<::reco::CALayersSoA>(),
0456 geometry_d.view<::reco::CAGraphSoA>(),
0457 queue);
0458
0459 HelixFit fitter(bfield, m_params.algoParams_.fitNas4_);
0460 fitter.allocate(kernels.tupleMultiplicity(), tracks.view(), kernels.hitContainer());
0461 if (m_params.algoParams_.useRiemannFit_) {
0462 fitter.launchRiemannKernels(hits_d.view(),
0463 geometry_d.view<::reco::CAModulesSoA>(),
0464 hits_d.view().metadata().size(),
0465 TrackerTraits::maxNumberOfQuadruplets,
0466 queue);
0467 } else {
0468 fitter.launchBrokenLineKernels(hits_d.view(),
0469 geometry_d.view<::reco::CAModulesSoA>(),
0470 hits_d.view().metadata().size(),
0471 TrackerTraits::maxNumberOfQuadruplets,
0472 queue);
0473 }
0474 kernels.classifyTuples(hits_d.view(), tracks.view(), queue);
0475 #ifdef GPU_DEBUG
0476 alpaka::wait(queue);
0477 std::cout << "finished building pixel tracks on GPU" << std::endl;
0478 #endif
0479
0480 return tracks;
0481 }
0482
0483 template class CAHitNtupletGenerator<pixelTopology::Phase1>;
0484 template class CAHitNtupletGenerator<pixelTopology::Phase2>;
0485 template class CAHitNtupletGenerator<pixelTopology::HIonPhase1>;
0486 }