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 <set>
5 #include <vector>
6 
7 #include <alpaka/alpaka.hpp>
8 
17 
19 
21  ALPAKA_FN_HOST_ACC ALPAKA_FN_ACC ALPAKA_FN_INLINE static bool hasQuality(const SiPixelMappingSoAConstView& view) {
22  return view.hasQuality();
23  }
24 
25  ALPAKA_FN_HOST_ACC ALPAKA_FN_ACC ALPAKA_FN_INLINE static cms::alpakatools::device_buffer<Device, unsigned char[]>
26  getModToUnpRegionalAsync(std::set<unsigned int> const& modules,
27  const SiPixelFedCablingTree* cabling,
28  std::vector<unsigned int> const& fedIds,
29  Queue& queue) {
30  auto modToUnpDevice = cms::alpakatools::make_device_buffer<unsigned char[]>(queue, pixelgpudetails::MAX_SIZE);
31  auto modToUnpHost = cms::alpakatools::make_host_buffer<unsigned char[]>(queue, pixelgpudetails::MAX_SIZE);
32 
33  unsigned int startFed = fedIds.front();
34  unsigned int endFed = fedIds.back() - 1;
35 
37  int index = 1;
38 
39  for (unsigned int fed = startFed; fed <= endFed; fed++) {
40  for (unsigned int link = 1; link <= pixelgpudetails::MAX_LINK; link++) {
41  for (unsigned int roc = 1; roc <= pixelgpudetails::MAX_ROC; roc++) {
42  path = {fed, link, roc};
43  const sipixelobjects::PixelROC* pixelRoc = cabling->findItem(path);
44  if (pixelRoc != nullptr) {
45  modToUnpHost[index] = (not modules.empty()) and (modules.find(pixelRoc->rawId()) == modules.end());
46  } else { // store some dummy number
47  modToUnpHost[index] = true;
48  }
49  index++;
50  }
51  }
52  }
53 
54  alpaka::memcpy(queue, modToUnpDevice, modToUnpHost);
55 
56  return modToUnpDevice;
57  }
58  };
59 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
60 #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:177
constexpr unsigned int MAX_LINK
uint32_t rawId() const
return the DetUnit to which this ROC belongs to.
Definition: PixelROC.h:34