Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-13 22:52:23

0001 /**
0002    Simple test for the reco::ZVertexSoA data structure
0003    which inherits from Portable{Host}Collection.
0004 
0005    Creates an instance of the class (automatically allocates
0006    memory on device), passes the view of the SoA data to
0007    the kernels which:
0008    - Fill the SoA with data.
0009    - Verify that the data written is correct.
0010 
0011    Then, the SoA data are copied back to Host, where
0012    a temporary host-side view (tmp_view) is created using
0013    the same Layout to access the data on host and print it.
0014  */
0015 
0016 #include <cstdlib>
0017 #include <unistd.h>
0018 
0019 #include <alpaka/alpaka.hpp>
0020 
0021 #include "DataFormats/VertexSoA/interface/ZVertexHost.h"
0022 #include "DataFormats/VertexSoA/interface/alpaka/ZVertexSoACollection.h"
0023 #include "FWCore/Utilities/interface/stringize.h"
0024 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0025 #include "HeterogeneousCore/AlpakaInterface/interface/devices.h"
0026 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0027 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0028 
0029 #include "ZVertexSoA_test.h"
0030 
0031 using namespace ALPAKA_ACCELERATOR_NAMESPACE;
0032 
0033 // Run 3 values, used for testing
0034 constexpr uint32_t maxTracks = 32 * 1024;
0035 constexpr uint32_t maxVertices = 1024;
0036 
0037 int main() {
0038   // Get the list of devices on the current platform
0039   auto const& devices = cms::alpakatools::devices<Platform>();
0040   if (devices.empty()) {
0041     std::cerr << "No devices available for the " EDM_STRINGIZE(ALPAKA_ACCELERATOR_NAMESPACE) " backend, "
0042       "the test will be skipped.\n";
0043     exit(EXIT_FAILURE);
0044   }
0045 
0046   // Run the test on each device
0047   for (const auto& device : devices) {
0048     Queue queue(device);
0049 
0050     // Inner scope to deallocate memory before destroying the stream
0051     {
0052       // Instantiate vertices on device. PortableCollection allocates
0053       // SoA on device automatically.
0054       ZVertexSoACollection zvertex_d({{maxTracks, maxVertices}}, queue);
0055       testZVertexSoAT::runKernels(zvertex_d.view(), zvertex_d.view<reco::ZVertexTracksSoA>(), queue);
0056 
0057       // If the device is actually the host, use the collection as-is.
0058       // Otherwise, copy the data from the device to the host.
0059       ZVertexHost zvertex_h;
0060 #ifdef ALPAKA_ACC_CPU_B_SEQ_T_SEQ_ENABLED
0061       zvertex_h = std::move(zvertex_d);
0062 #else
0063       zvertex_h = cms::alpakatools::CopyToHost<ZVertexSoACollection>::copyAsync(queue, zvertex_d);
0064 #endif
0065       alpaka::wait(queue);
0066       std::cout << zvertex_h.view().metadata().size() << std::endl;
0067 
0068       // Print results
0069       std::cout << "idv\t"
0070                 << "zv\t"
0071                 << "wv\t"
0072                 << "chi2\t"
0073                 << "ptv2\t"
0074                 << "ndof\t"
0075                 << "sortInd\t"
0076                 << "nvFinal\n";
0077 
0078       auto vtx_v = zvertex_h.view<reco::ZVertexSoA>();
0079       auto trk_v = zvertex_h.view<reco::ZVertexTracksSoA>();
0080       for (int i = 0; i < 10; ++i) {
0081         auto vi = vtx_v[i];
0082         auto ti = trk_v[i];
0083         std::cout << (int)ti.idv() << "\t" << vi.zv() << "\t" << vi.wv() << "\t" << vi.chi2() << "\t" << vi.ptv2()
0084                   << "\t" << (int)ti.ndof() << "\t" << vi.sortInd() << "\t" << (int)vtx_v.nvFinal() << std::endl;
0085       }
0086     }
0087   }
0088 
0089   return EXIT_SUCCESS;
0090 }