Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:40

0001 #ifndef RecoTracker_PixelVertexFinding_plugins_alpaka_vertexFinder_h
0002 #define RecoTracker_PixelVertexFinding_plugins_alpaka_vertexFinder_h
0003 
0004 #include <cstddef>
0005 #include <cstdint>
0006 
0007 #include <alpaka/alpaka.hpp>
0008 
0009 #include "DataFormats/TrackSoA/interface/alpaka/TrackUtilities.h"
0010 #include "DataFormats/VertexSoA/interface/ZVertexDevice.h"
0011 #include "DataFormats/VertexSoA/interface/ZVertexHost.h"
0012 #include "DataFormats/VertexSoA/interface/ZVertexSoA.h"
0013 #include "DataFormats/VertexSoA/interface/alpaka/ZVertexSoACollection.h"
0014 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0015 #include "RecoTracker/PixelVertexFinding/interface/PixelVertexWorkSpaceLayout.h"
0016 #include "RecoTracker/PixelVertexFinding/plugins/alpaka/PixelVertexWorkSpaceSoADeviceAlpaka.h"
0017 
0018 namespace ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder {
0019 
0020   using namespace cms::alpakatools;
0021   using VtxSoAView = ::reco::ZVertexSoAView;
0022   using WsSoAView = ::vertexFinder::PixelVertexWorkSpaceSoAView;
0023 
0024   class Init {
0025   public:
0026     template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
0027     ALPAKA_FN_ACC void operator()(const TAcc &acc, VtxSoAView pdata, WsSoAView pws) const {
0028       pdata.nvFinal() = 0;  // initialization
0029       ::vertexFinder::init(pws);
0030     }
0031   };
0032 
0033   template <typename TrackerTraits>
0034   class Producer {
0035     using TkSoAConstView = reco::TrackSoAConstView<TrackerTraits>;
0036 
0037   public:
0038     Producer(bool oneKernel,
0039              bool useDensity,
0040              bool useDBSCAN,
0041              bool useIterative,
0042              bool doSplitting,
0043              int iminT,      // min number of neighbours to be "core"
0044              float ieps,     // max absolute distance to cluster
0045              float ierrmax,  // max error to be "seed"
0046              float ichi2max  // max normalized distance to cluster
0047              )
0048         : oneKernel_(oneKernel && !(useDBSCAN || useIterative)),
0049           useDensity_(useDensity),
0050           useDBSCAN_(useDBSCAN),
0051           useIterative_(useIterative),
0052           doSplitting_(doSplitting),
0053           minT(iminT),
0054           eps(ieps),
0055           errmax(ierrmax),
0056           chi2max(ichi2max) {}
0057 
0058     ~Producer() = default;
0059 
0060     ZVertexSoACollection makeAsync(Queue &queue, const TkSoAConstView &tracks_view, float ptMin, float ptMax) const;
0061 
0062   private:
0063     const bool oneKernel_;     // run everything (cluster,fit,split,sort) in one kernel. Uses only density clusterizer
0064     const bool useDensity_;    // use density clusterizer
0065     const bool useDBSCAN_;     // use DBScan clusterizer
0066     const bool useIterative_;  // use iterative clusterizer
0067     const bool doSplitting_;   //run vertex splitting
0068 
0069     int minT;       // min number of neighbours to be "core"
0070     float eps;      // max absolute distance to cluster
0071     float errmax;   // max error to be "seed"
0072     float chi2max;  // max normalized distance to cluster
0073   };
0074 
0075 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder
0076 
0077 #endif  // RecoTracker_PixelVertexFinding_plugins_alpaka_vertexFinder_h