CMS 3D CMS Logo

List of all members | Public Member Functions
ALPAKA_ACCELERATOR_NAMESPACE::FastCluster Class Reference

Public Member Functions

template<bool debug = false, typename TAcc , typename = std::enable_if<!std::is_same_v<Device, alpaka::DevCpu>>>
ALPAKA_FN_ACC void operator() (const TAcc &acc, const reco::PFRecHitHostCollection::ConstView pfRecHits, const reco::PFClusterParamsDeviceCollection::ConstView pfClusParams, const reco::PFRecHitHCALTopologyDeviceCollection::ConstView topology, reco::PFClusteringVarsDeviceCollection::View pfClusteringVars, reco::PFClusterDeviceCollection::View clusterView, reco::PFRecHitFractionDeviceCollection::View fracView) const
 

Detailed Description

Definition at line 1349 of file PFClusterSoAProducerKernel.dev.cc.

Member Function Documentation

◆ operator()()

template<bool debug = false, typename TAcc , typename = std::enable_if<!std::is_same_v<Device, alpaka::DevCpu>>>
ALPAKA_FN_ACC void ALPAKA_ACCELERATOR_NAMESPACE::FastCluster::operator() ( const TAcc &  acc,
const reco::PFRecHitHostCollection::ConstView  pfRecHits,
const reco::PFClusterParamsDeviceCollection::ConstView  pfClusParams,
const reco::PFRecHitHCALTopologyDeviceCollection::ConstView  topology,
reco::PFClusteringVarsDeviceCollection::View  pfClusteringVars,
reco::PFClusterDeviceCollection::View  clusterView,
reco::PFRecHitFractionDeviceCollection::View  fracView 
) const
inline

Definition at line 1352 of file PFClusterSoAProducerKernel.dev.cc.

References ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), debug, ALPAKA_ACCELERATOR_NAMESPACE::hcalFastCluster_multiSeedIterative(), ALPAKA_ACCELERATOR_NAMESPACE::hcalFastCluster_multiSeedParallel(), ALPAKA_ACCELERATOR_NAMESPACE::hcalFastCluster_singleSeed(), cms::alpakatools::once_per_block(), HLT_2024v14_cff::pfRecHits, ALPAKA_ACCELERATOR_NAMESPACE::threadsPerBlockForClustering, and HLT_2024v14_cff::topology.

1358  {
1359  const int nRH = pfRecHits.size();
1360  int& topoId = alpaka::declareSharedVar<int, __COUNTER__>(acc);
1361  int& nRHTopo = alpaka::declareSharedVar<int, __COUNTER__>(acc);
1362  int& nSeeds = alpaka::declareSharedVar<int, __COUNTER__>(acc);
1363 
1364  if (once_per_block(acc)) {
1365  topoId = alpaka::getIdx<alpaka::Grid, alpaka::Blocks>(acc)[0u];
1366  nRHTopo = pfClusteringVars[topoId].topoRHCount();
1367  nSeeds = pfClusteringVars[topoId].topoSeedCount();
1368  }
1369 
1370  alpaka::syncBlockThreads(acc); // all threads call sync
1371 
1372  if (topoId < nRH && nRHTopo > 0 && nSeeds > 0) {
1373  if (nRHTopo == nSeeds) {
1374  // PF cluster is isolated seed. No iterations needed
1375  if (once_per_block(acc)) {
1376  // Fill PFCluster-level information
1377  int rhIdx = pfClusteringVars[pfClusteringVars[topoId].topoSeedOffsets()].topoSeedList();
1378  int seedIdx = pfClusteringVars[rhIdx].rhIdxToSeedIdx();
1379  clusterView[seedIdx].energy() = pfRecHits[rhIdx].energy();
1380  clusterView[seedIdx].x() = pfRecHits[rhIdx].x();
1381  clusterView[seedIdx].y() = pfRecHits[rhIdx].y();
1382  clusterView[seedIdx].z() = pfRecHits[rhIdx].z();
1383  }
1384  // singleSeed and multiSeedParallel functions work only for GPU backend
1385  } else if ((not std::is_same_v<Device, alpaka::DevCpu>)&&nSeeds == 1) {
1386  // Single seed cluster
1388  acc, pfClusParams, topology, topoId, nRHTopo, pfRecHits, pfClusteringVars, clusterView, fracView);
1389  } else if ((not std::is_same_v<Device, alpaka::DevCpu>)&&nSeeds <= 100 &&
1390  nRHTopo - nSeeds < threadsPerBlockForClustering) {
1392  acc, pfClusParams, topology, topoId, nSeeds, nRHTopo, pfRecHits, pfClusteringVars, clusterView, fracView);
1393  } else if (nSeeds <= 400 && nRHTopo - nSeeds <= 1500) {
1394  // nSeeds value must match exotic in FastClusterExotic
1396  acc, pfClusParams, topology, topoId, nSeeds, nRHTopo, pfRecHits, pfClusteringVars, clusterView, fracView);
1397  } else {
1398  if constexpr (debug) {
1399  if (once_per_block(acc))
1400  printf("Topo cluster %d has %d seeds and %d rechits. Will be processed in next kernel.\n",
1401  topoId,
1402  nSeeds,
1403  nRHTopo);
1404  }
1405  }
1406  }
1407  }
ALPAKA_FN_ACC constexpr bool once_per_block(TAcc const &acc)
static ALPAKA_FN_ACC void hcalFastCluster_multiSeedIterative(const TAcc &acc, reco::PFClusterParamsDeviceCollection::ConstView pfClusParams, const reco::PFRecHitHCALTopologyDeviceCollection::ConstView topology, int topoId, int nSeeds, int nRHTopo, reco::PFRecHitDeviceCollection::ConstView pfRecHits, reco::PFClusteringVarsDeviceCollection::View pfClusteringVars, reco::PFClusterDeviceCollection::View clusterView, reco::PFRecHitFractionDeviceCollection::View fracView)
static ALPAKA_FN_ACC void hcalFastCluster_multiSeedParallel(const TAcc &acc, reco::PFClusterParamsDeviceCollection::ConstView pfClusParams, const reco::PFRecHitHCALTopologyDeviceCollection::ConstView topology, int topoId, int nSeeds, int nRHTopo, reco::PFRecHitDeviceCollection::ConstView pfRecHits, reco::PFClusteringVarsDeviceCollection::View pfClusteringVars, reco::PFClusterDeviceCollection::View clusterView, reco::PFRecHitFractionDeviceCollection::View fracView)
#define debug
Definition: HDRShower.cc:19
static ALPAKA_FN_ACC void hcalFastCluster_singleSeed(const TAcc &acc, reco::PFClusterParamsDeviceCollection::ConstView pfClusParams, const reco::PFRecHitHCALTopologyDeviceCollection::ConstView topology, int topoId, int nRHTopo, reco::PFRecHitDeviceCollection::ConstView pfRecHits, reco::PFClusteringVarsDeviceCollection::View pfClusteringVars, reco::PFClusterDeviceCollection::View clusterView, reco::PFRecHitFractionDeviceCollection::View fracView)