Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-01-19 02:53:21

0001 #ifndef RecoPixelVertexing_PixelTriplets_plugins_CAHitNtupletGeneratorKernels_h
0002 #define RecoPixelVertexing_PixelTriplets_plugins_CAHitNtupletGeneratorKernels_h
0003 
0004 // #define GPU_DEBUG
0005 
0006 #include "GPUCACell.h"
0007 #include "gpuPixelDoublets.h"
0008 
0009 #include "CUDADataFormats/Track/interface/PixelTrackUtilities.h"
0010 #include "CUDADataFormats/TrackingRecHit/interface/TrackingRecHitsUtilities.h"
0011 #include "CUDADataFormats/Common/interface/HeterogeneousSoA.h"
0012 #include "CUDADataFormats/TrackingRecHit/interface/TrackingRecHitSoADevice.h"
0013 #include "CUDADataFormats/Track/interface/TrackSoAHeterogeneousHost.h"
0014 
0015 // #define DUMP_GPU_TK_TUPLES
0016 
0017 namespace caHitNtupletGenerator {
0018 
0019   //Configuration params common to all topologies, for the algorithms
0020   struct AlgoParams {
0021     const bool onGPU_;
0022     const uint32_t minHitsForSharingCut_;
0023     const bool useRiemannFit_;
0024     const bool fitNas4_;
0025     const bool includeJumpingForwardDoublets_;
0026     const bool earlyFishbone_;
0027     const bool lateFishbone_;
0028     const bool doStats_;
0029     const bool doSharedHitCut_;
0030     const bool dupPassThrough_;
0031     const bool useSimpleTripletCleaner_;
0032   };
0033 
0034   //CAParams
0035   struct CACommon {
0036     const uint32_t minHitsPerNtuplet_;
0037     const float ptmin_;
0038     const float CAThetaCutBarrel_;
0039     const float CAThetaCutForward_;
0040     const float hardCurvCut_;
0041     const float dcaCutInnerTriplet_;
0042     const float dcaCutOuterTriplet_;
0043   };
0044 
0045   template <typename TrackerTraits, typename Enable = void>
0046   struct CAParamsT : public CACommon {
0047     __device__ __forceinline__ bool startingLayerPair(int16_t pid) const { return false; };
0048     __device__ __forceinline__ bool startAt0(int16_t pid) const { return false; };
0049   };
0050 
0051   template <typename TrackerTraits>
0052   struct CAParamsT<TrackerTraits, pixelTopology::isPhase1Topology<TrackerTraits>> : public CACommon {
0053     /// Is is a starting layer pair?
0054     __device__ __forceinline__ bool startingLayerPair(int16_t pid) const {
0055       return minHitsPerNtuplet_ > 3 ? pid < 3 : pid < 8 || pid > 12;
0056     }
0057 
0058     /// Is this a pair with inner == 0?
0059     __device__ __forceinline__ bool startAt0(int16_t pid) const {
0060       assert((pixelTopology::Phase1::layerPairs[pid * 2] == 0) ==
0061              (pid < 3 || pid == 13 || pid == 15 || pid == 16));  // to be 100% sure it's working, may be removed
0062       return pixelTopology::Phase1::layerPairs[pid * 2] == 0;
0063     }
0064   };
0065 
0066   template <typename TrackerTraits>
0067   struct CAParamsT<TrackerTraits, pixelTopology::isPhase2Topology<TrackerTraits>> : public CACommon {
0068     const bool includeFarForwards_;
0069     /// Is is a starting layer pair?
0070     __device__ __forceinline__ bool startingLayerPair(int16_t pid) const {
0071       return pid < 33;  // in principle one could remove 5,6,7 23, 28 and 29
0072     }
0073 
0074     /// Is this a pair with inner == 0
0075     __device__ __forceinline__ bool startAt0(int16_t pid) const {
0076       assert((pixelTopology::Phase2::layerPairs[pid * 2] == 0) == ((pid < 3) | (pid >= 23 && pid < 28)));
0077       return pixelTopology::Phase2::layerPairs[pid * 2] == 0;
0078     }
0079   };
0080 
0081   //Full list of params = algo params + ca params + cell params + quality cuts
0082   //Generic template
0083   template <typename TrackerTraits, typename Enable = void>
0084   struct ParamsT : public AlgoParams {
0085     // one should define the params for its own pixelTopology
0086     // not defining anything here
0087     inline uint32_t nPairs() const { return 0; }
0088   };
0089 
0090   template <typename TrackerTraits>
0091   struct ParamsT<TrackerTraits, pixelTopology::isPhase1Topology<TrackerTraits>> : public AlgoParams {
0092     using TT = TrackerTraits;
0093     using QualityCuts = pixelTrack::QualityCutsT<TT>;  //track quality cuts
0094     using CellCuts = gpuPixelDoublets::CellCutsT<TT>;  //cell building cuts
0095     using CAParams = CAParamsT<TT>;                    //params to be used on device
0096 
0097     ParamsT(AlgoParams const& commonCuts,
0098             CellCuts const& cellCuts,
0099             QualityCuts const& cutsCuts,
0100             CAParams const& caParams)
0101         : AlgoParams(commonCuts), cellCuts_(cellCuts), qualityCuts_(cutsCuts), caParams_(caParams) {}
0102 
0103     const CellCuts cellCuts_;
0104     const QualityCuts qualityCuts_{// polynomial coefficients for the pT-dependent chi2 cut
0105                                    {0.68177776, 0.74609577, -0.08035491, 0.00315399},
0106                                    // max pT used to determine the chi2 cut
0107                                    10.,
0108                                    // chi2 scale factor: 30 for broken line fit, 45 for Riemann fit
0109                                    30.,
0110                                    // regional cuts for triplets
0111                                    {
0112                                        0.3,  // |Tip| < 0.3 cm
0113                                        0.5,  // pT > 0.5 GeV
0114                                        12.0  // |Zip| < 12.0 cm
0115                                    },
0116                                    // regional cuts for quadruplets
0117                                    {
0118                                        0.5,  // |Tip| < 0.5 cm
0119                                        0.3,  // pT > 0.3 GeV
0120                                        12.0  // |Zip| < 12.0 cm
0121                                    }};
0122     const CAParams caParams_;
0123     /// Compute the number of pairs
0124     inline uint32_t nPairs() const {
0125       // take all layer pairs into account
0126       uint32_t nActualPairs = TT::nPairs;
0127       if (not includeJumpingForwardDoublets_) {
0128         // exclude forward "jumping" layer pairs
0129         nActualPairs = TT::nPairsForTriplets;
0130       }
0131       if (caParams_.minHitsPerNtuplet_ > 3) {
0132         // for quadruplets, exclude all "jumping" layer pairs
0133         nActualPairs = TT::nPairsForQuadruplets;
0134       }
0135 
0136       return nActualPairs;
0137     }
0138 
0139   };  // Params Phase1
0140 
0141   template <typename TrackerTraits>
0142   struct ParamsT<TrackerTraits, pixelTopology::isPhase2Topology<TrackerTraits>> : public AlgoParams {
0143     using TT = TrackerTraits;
0144     using QualityCuts = pixelTrack::QualityCutsT<TT>;
0145     using CellCuts = gpuPixelDoublets::CellCutsT<TT>;
0146     using CAParams = CAParamsT<TT>;
0147 
0148     ParamsT(AlgoParams const& commonCuts,
0149             CellCuts const& cellCuts,
0150             QualityCuts const& qualityCuts,
0151             CAParams const& caParams)
0152         : AlgoParams(commonCuts), cellCuts_(cellCuts), qualityCuts_(qualityCuts), caParams_(caParams) {}
0153 
0154     // quality cuts
0155     const CellCuts cellCuts_;
0156     const QualityCuts qualityCuts_{5.0f, /*chi2*/ 0.9f, /* pT in Gev*/ 0.4f, /*zip in cm*/ 12.0f /*tip in cm*/};
0157     const CAParams caParams_;
0158 
0159     inline uint32_t nPairs() const {
0160       // take all layer pairs into account
0161       uint32_t nActualPairs = TT::nPairsMinimal;
0162       if (caParams_.includeFarForwards_) {
0163         // considera far forwards (> 11 & > 23)
0164         nActualPairs = TT::nPairsFarForwards;
0165       }
0166       if (includeJumpingForwardDoublets_) {
0167         // include jumping forwards
0168         nActualPairs = TT::nPairs;
0169       }
0170 
0171       return nActualPairs;
0172     }
0173 
0174   };  // Params Phase1
0175 
0176   // counters
0177   struct Counters {
0178     unsigned long long nEvents;
0179     unsigned long long nHits;
0180     unsigned long long nCells;
0181     unsigned long long nTuples;
0182     unsigned long long nFitTracks;
0183     unsigned long long nLooseTracks;
0184     unsigned long long nGoodTracks;
0185     unsigned long long nUsedHits;
0186     unsigned long long nDupHits;
0187     unsigned long long nFishCells;
0188     unsigned long long nKilledCells;
0189     unsigned long long nEmptyCells;
0190     unsigned long long nZeroTrackCells;
0191   };
0192 
0193   using Quality = pixelTrack::Quality;
0194 
0195 }  // namespace caHitNtupletGenerator
0196 
0197 template <typename TTraits, typename TTTraits>
0198 class CAHitNtupletGeneratorKernels {
0199 public:
0200   using Traits = TTraits;
0201   using TrackerTraits = TTTraits;
0202   using QualityCuts = pixelTrack::QualityCutsT<TrackerTraits>;
0203   using Params = caHitNtupletGenerator::ParamsT<TrackerTraits>;
0204   using CAParams = caHitNtupletGenerator::CAParamsT<TrackerTraits>;
0205   using Counters = caHitNtupletGenerator::Counters;
0206 
0207   template <typename T>
0208   using unique_ptr = typename Traits::template unique_ptr<T>;
0209 
0210   using HitsView = TrackingRecHitSoAView<TrackerTraits>;
0211   using HitsConstView = TrackingRecHitSoAConstView<TrackerTraits>;
0212   using TkSoAView = TrackSoAView<TrackerTraits>;
0213 
0214   using HitToTuple = caStructures::HitToTupleT<TrackerTraits>;
0215   using TupleMultiplicity = caStructures::TupleMultiplicityT<TrackerTraits>;
0216   using CellNeighborsVector = caStructures::CellNeighborsVectorT<TrackerTraits>;
0217   using CellNeighbors = caStructures::CellNeighborsT<TrackerTraits>;
0218   using CellTracksVector = caStructures::CellTracksVectorT<TrackerTraits>;
0219   using CellTracks = caStructures::CellTracksT<TrackerTraits>;
0220   using OuterHitOfCellContainer = caStructures::OuterHitOfCellContainerT<TrackerTraits>;
0221   using OuterHitOfCell = caStructures::OuterHitOfCellT<TrackerTraits>;
0222 
0223   using CACell = GPUCACellT<TrackerTraits>;
0224 
0225   using Quality = pixelTrack::Quality;
0226   using HitContainer = typename TrackSoA<TrackerTraits>::HitContainer;
0227 
0228   CAHitNtupletGeneratorKernels(Params const& params)
0229       : params_(params), paramsMaxDoubletes3Quarters_(3 * params.cellCuts_.maxNumberOfDoublets_ / 4) {}
0230 
0231   ~CAHitNtupletGeneratorKernels() = default;
0232 
0233   TupleMultiplicity const* tupleMultiplicity() const { return device_tupleMultiplicity_.get(); }
0234 
0235   void launchKernels(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
0236 
0237   void classifyTuples(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
0238 
0239   void buildDoublets(const HitsConstView& hh, int32_t offsetBPIX2, cudaStream_t stream);
0240   void allocateOnGPU(int32_t nHits, cudaStream_t stream);
0241   void cleanup(cudaStream_t cudaStream);
0242 
0243   static void printCounters(Counters const* counters);
0244   void setCounters(Counters* counters) { counters_ = counters; }
0245 
0246 protected:
0247   Counters* counters_ = nullptr;
0248   // workspace
0249   unique_ptr<unsigned char[]> cellStorage_;
0250   unique_ptr<CellNeighborsVector> device_theCellNeighbors_;
0251   CellNeighbors* device_theCellNeighborsContainer_;
0252   unique_ptr<CellTracksVector> device_theCellTracks_;
0253   CellTracks* device_theCellTracksContainer_;
0254 
0255   unique_ptr<CACell[]> device_theCells_;
0256   unique_ptr<OuterHitOfCellContainer[]> device_isOuterHitOfCell_;
0257   OuterHitOfCell isOuterHitOfCell_;
0258   uint32_t* device_nCells_ = nullptr;
0259 
0260   unique_ptr<HitToTuple> device_hitToTuple_;
0261   unique_ptr<uint32_t[]> device_hitToTupleStorage_;
0262   typename HitToTuple::View hitToTupleView_;
0263 
0264   cms::cuda::AtomicPairCounter* device_hitToTuple_apc_ = nullptr;
0265 
0266   cms::cuda::AtomicPairCounter* device_hitTuple_apc_ = nullptr;
0267 
0268   unique_ptr<TupleMultiplicity> device_tupleMultiplicity_;
0269 
0270   unique_ptr<cms::cuda::AtomicPairCounter::c_type[]> device_storage_;
0271 
0272   // params
0273   Params params_;
0274   /// Intermediate result avoiding repeated computations.
0275   const uint32_t paramsMaxDoubletes3Quarters_;
0276   /// Compute the number of doublet blocks for block size
0277   inline uint32_t nDoubletBlocks(uint32_t blockSize) {
0278     // We want (3 * params_.maxNumberOfDoublets_ / 4 + blockSize - 1) / blockSize, but first part is pre-computed.
0279     return (paramsMaxDoubletes3Quarters_ + blockSize - 1) / blockSize;
0280   }
0281 
0282   /// Compute the number of quadruplet blocks for block size
0283   inline uint32_t nQuadrupletBlocks(uint32_t blockSize) {
0284     // pixelTopology::maxNumberOfQuadruplets is a constexpr, so the compiler will pre compute the 3*max/4
0285     return (3 * TrackerTraits::maxNumberOfQuadruplets / 4 + blockSize - 1) / blockSize;
0286   }
0287 };
0288 
0289 template <typename TrackerTraits>
0290 class CAHitNtupletGeneratorKernelsGPU : public CAHitNtupletGeneratorKernels<cms::cudacompat::GPUTraits, TrackerTraits> {
0291   using CAHitNtupletGeneratorKernels<cms::cudacompat::GPUTraits, TrackerTraits>::CAHitNtupletGeneratorKernels;
0292 
0293   using Counters = caHitNtupletGenerator::Counters;
0294   using CAParams = caHitNtupletGenerator::CAParamsT<TrackerTraits>;
0295 
0296   using HitContainer = typename TrackSoA<TrackerTraits>::HitContainer;
0297 
0298   using CellNeighborsVector = caStructures::CellNeighborsVectorT<TrackerTraits>;
0299   using HitToTuple = caStructures::HitToTupleT<TrackerTraits>;
0300   using CellTracksVector = caStructures::CellTracksVectorT<TrackerTraits>;
0301   using TupleMultiplicity = caStructures::TupleMultiplicityT<TrackerTraits>;
0302 
0303   using HitsConstView = TrackingRecHitSoAConstView<TrackerTraits>;
0304   using TkSoAView = TrackSoAView<TrackerTraits>;
0305 
0306 public:
0307   void launchKernels(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
0308   void classifyTuples(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
0309   void buildDoublets(const HitsConstView& hh, int32_t offsetBPIX2, cudaStream_t stream);
0310   void allocateOnGPU(int32_t nHits, cudaStream_t stream);
0311   static void printCounters(Counters const* counters);
0312 };
0313 
0314 template <typename TrackerTraits>
0315 class CAHitNtupletGeneratorKernelsCPU : public CAHitNtupletGeneratorKernels<cms::cudacompat::CPUTraits, TrackerTraits> {
0316   using CAHitNtupletGeneratorKernels<cms::cudacompat::CPUTraits, TrackerTraits>::CAHitNtupletGeneratorKernels;
0317 
0318   using Counters = caHitNtupletGenerator::Counters;
0319   using CAParams = caHitNtupletGenerator::CAParamsT<TrackerTraits>;
0320 
0321   using HitContainer = typename TrackSoA<TrackerTraits>::HitContainer;
0322 
0323   using CellNeighborsVector = caStructures::CellNeighborsVectorT<TrackerTraits>;
0324   using HitToTuple = caStructures::HitToTupleT<TrackerTraits>;
0325   using CellTracksVector = caStructures::CellTracksVectorT<TrackerTraits>;
0326   using TupleMultiplicity = caStructures::TupleMultiplicityT<TrackerTraits>;
0327 
0328   using HitsConstView = TrackingRecHitSoAConstView<TrackerTraits>;
0329   using TkSoAView = TrackSoAView<TrackerTraits>;
0330 
0331 public:
0332   void launchKernels(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
0333   void classifyTuples(const HitsConstView& hh, TkSoAView& track_view, cudaStream_t cudaStream);
0334   void buildDoublets(const HitsConstView& hh, int32_t offsetBPIX2, cudaStream_t stream);
0335   void allocateOnGPU(int32_t nHits, cudaStream_t stream);
0336   static void printCounters(Counters const* counters);
0337 };
0338 
0339 #endif  // RecoPixelVertexing_PixelTriplets_plugins_CAHitNtupletGeneratorKernels_h