File indexing completed on 2024-04-06 12:02:33
0001 #ifndef CondFormats_SiPixelObjects_interface_alpaka_SiPixelMappingUtilities_h
0002 #define CondFormats_SiPixelObjects_interface_alpaka_SiPixelMappingUtilities_h
0003
0004 #include <set>
0005 #include <vector>
0006
0007 #include <alpaka/alpaka.hpp>
0008
0009 #include "CondFormats/SiPixelObjects/interface/CablingPathToDetUnit.h"
0010 #include "CondFormats/SiPixelObjects/interface/PixelROC.h"
0011 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
0012 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingTree.h"
0013 #include "CondFormats/SiPixelObjects/interface/SiPixelMappingLayout.h"
0014 #include "CondFormats/SiPixelObjects/interface/SiPixelROCsStatusAndMapping.h"
0015 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0016 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0017
0018 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0019
0020 struct SiPixelMappingUtilities {
0021 ALPAKA_FN_HOST_ACC ALPAKA_FN_ACC ALPAKA_FN_INLINE static bool hasQuality(const SiPixelMappingSoAConstView& view) {
0022 return view.hasQuality();
0023 }
0024
0025 ALPAKA_FN_HOST_ACC ALPAKA_FN_ACC ALPAKA_FN_INLINE static cms::alpakatools::device_buffer<Device, unsigned char[]>
0026 getModToUnpRegionalAsync(std::set<unsigned int> const& modules,
0027 const SiPixelFedCablingTree* cabling,
0028 std::vector<unsigned int> const& fedIds,
0029 Queue& queue) {
0030 auto modToUnpDevice = cms::alpakatools::make_device_buffer<unsigned char[]>(queue, pixelgpudetails::MAX_SIZE);
0031 auto modToUnpHost = cms::alpakatools::make_host_buffer<unsigned char[]>(queue, pixelgpudetails::MAX_SIZE);
0032
0033 unsigned int startFed = fedIds.front();
0034 unsigned int endFed = fedIds.back() - 1;
0035
0036 sipixelobjects::CablingPathToDetUnit path;
0037 int index = 1;
0038
0039 for (unsigned int fed = startFed; fed <= endFed; fed++) {
0040 for (unsigned int link = 1; link <= pixelgpudetails::MAX_LINK; link++) {
0041 for (unsigned int roc = 1; roc <= pixelgpudetails::MAX_ROC; roc++) {
0042 path = {fed, link, roc};
0043 const sipixelobjects::PixelROC* pixelRoc = cabling->findItem(path);
0044 if (pixelRoc != nullptr) {
0045 modToUnpHost[index] = (not modules.empty()) and (modules.find(pixelRoc->rawId()) == modules.end());
0046 } else {
0047 modToUnpHost[index] = true;
0048 }
0049 index++;
0050 }
0051 }
0052 }
0053
0054 alpaka::memcpy(queue, modToUnpDevice, modToUnpHost);
0055
0056 return modToUnpDevice;
0057 }
0058 };
0059 }
0060 #endif