CMS 3D CMS Logo

SiPixelMappingUtilities.h
Go to the documentation of this file.
1 #ifndef CondFormats_SiPixelObjects_interface_alpaka_SiPixelMappingUtilities_h
2 #define CondFormats_SiPixelObjects_interface_alpaka_SiPixelMappingUtilities_h
3 
4 #include <cstdint>
5 #include <alpaka/alpaka.hpp>
10 
12 
14  ALPAKA_FN_HOST_ACC ALPAKA_FN_ACC ALPAKA_FN_INLINE static bool hasQuality(const SiPixelMappingSoAConstView& view) {
15  return view.hasQuality();
16  }
17 
18  ALPAKA_FN_HOST_ACC ALPAKA_FN_ACC ALPAKA_FN_INLINE static cms::alpakatools::device_buffer<Device, unsigned char[]>
19  getModToUnpRegionalAsync(std::set<unsigned int> const& modules,
20  const SiPixelFedCablingTree* cabling,
21  std::vector<unsigned int> const& fedIds,
22  Queue& queue) {
23  auto modToUnpDevice = cms::alpakatools::make_device_buffer<unsigned char[]>(queue, pixelgpudetails::MAX_SIZE);
24  auto modToUnpHost = cms::alpakatools::make_host_buffer<unsigned char[]>(queue, pixelgpudetails::MAX_SIZE);
25 
26  unsigned int startFed = fedIds.front();
27  unsigned int endFed = fedIds.back() - 1;
28 
30  int index = 1;
31 
32  for (unsigned int fed = startFed; fed <= endFed; fed++) {
33  for (unsigned int link = 1; link <= pixelgpudetails::MAX_LINK; link++) {
34  for (unsigned int roc = 1; roc <= pixelgpudetails::MAX_ROC; roc++) {
35  path = {fed, link, roc};
36  const sipixelobjects::PixelROC* pixelRoc = cabling->findItem(path);
37  if (pixelRoc != nullptr) {
38  modToUnpHost[index] = (not modules.empty()) and (modules.find(pixelRoc->rawId()) == modules.end());
39  } else { // store some dummy number
40  modToUnpHost[index] = true;
41  }
42  index++;
43  }
44  }
45  }
46 
47  alpaka::memcpy(queue, modToUnpDevice, modToUnpHost);
48 
49  return modToUnpDevice;
50  }
51  };
52 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
53 #endif //CondFormats_SiPixelObjects_interface_alpaka_SiPixelMappingUtilities_h
const sipixelobjects::PixelROC * findItem(const sipixelobjects::CablingPathToDetUnit &path) const final
constexpr unsigned int MAX_ROC
ALPAKA_FN_HOST_ACC ALPAKA_FN_ACC static ALPAKA_FN_INLINE bool hasQuality(const SiPixelMappingSoAConstView &view)
constexpr unsigned int MAX_SIZE
ALPAKA_FN_HOST_ACC ALPAKA_FN_ACC static ALPAKA_FN_INLINE cms::alpakatools::device_buffer< Device, unsigned char[]> getModToUnpRegionalAsync(std::set< unsigned int > const &modules, const SiPixelFedCablingTree *cabling, std::vector< unsigned int > const &fedIds, Queue &queue)
typename detail::buffer_type< TDev, T >::type device_buffer
Definition: memory.h:169
constexpr unsigned int MAX_LINK
uint32_t rawId() const
return the DetUnit to which this ROC belongs to.
Definition: PixelROC.h:34