Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-21 01:39:53

0001 #define CATCH_CONFIG_MAIN
0002 #include <catch.hpp>
0003 #include <numbers>
0004 #include <vector>
0005 
0006 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0007 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0008 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0009 #include "HeterogeneousCore/AlpakaMath/interface/deltaPhi.h"
0010 
0011 using namespace cms::alpakatools;
0012 using namespace ALPAKA_ACCELERATOR_NAMESPACE;
0013 
0014 struct phiFuncsUnitTestsKernel {
0015   template <typename TAcc, typename T>
0016   ALPAKA_FN_ACC void operator()(TAcc const& acc, T* out) const {
0017     // Unit circle typical values
0018     out[0] = phi<TAcc, T>(acc, 1.0, 0.0);          // x = 1.0, y = 0.0 => phi = 0
0019     out[1] = phi<TAcc, T>(acc, 0.0, 1.0);          // x = 0.0, y = 1.0 => phi = pi/2
0020     out[2] = phi<TAcc, T>(acc, -1.0, 0.0);         // x = -1.0, y = 0.0 => phi = pi
0021     out[3] = phi<TAcc, T>(acc, 0.0, -1.0);         // x = 0.0, y = -1.0 => phi = -pi/2
0022     out[4] = phi<TAcc, T>(acc, 0.7071, 0.7071);    // x = sqrt(2)/2, y = sqrt(2)/2 => phi = pi/4
0023     out[5] = phi<TAcc, T>(acc, -0.7071, -0.7071);  // x = sqrt(2)/2, y = sqrt(2)/2 => phi = -3pi/4
0024 
0025     // Making sure that delta phi is within [-pi, pi] range
0026     // Phi from unit circle
0027     out[6] = deltaPhi<TAcc, T>(acc, 1.0, 0.0, 0.0, -1.0);                // 3pi/2 - 0 = -pi/2
0028     out[7] = deltaPhi<TAcc, T>(acc, 0.0, 1.0, 0.0, -1.0);                // 3pi/2 - pi/2 = pi
0029     out[8] = deltaPhi<TAcc, T>(acc, 0.0, -1.0, 0.0, 1.0);                // pi/2 - 3pi/2 = -pi
0030     out[9] = deltaPhi<TAcc, T>(acc, 0.7071, -0.7071, -0.7071, -0.7071);  // -3pi/4 - (-pi/4) = -pi/2
0031 
0032     // Calculation directly from phi
0033     out[10] = deltaPhi<TAcc, T>(acc, 3. * M_PI / 2., 0.);           // 3pi/2 - 0 = -pi/2
0034     out[11] = deltaPhi<TAcc, T>(acc, 3. * M_PI / 2., M_PI / 2.);    // 3pi/2 - pi/2 = pi
0035     out[12] = deltaPhi<TAcc, T>(acc, M_PI / 2., 3. * M_PI / 2.);    // pi/2 - 3pi/2 = -pi
0036     out[13] = deltaPhi<TAcc, T>(acc, -3. * M_PI / 4., -M_PI / 4.);  // -3pi/4 - (-pi/4) = -pi/2
0037   }
0038 };
0039 
0040 template <typename T>
0041 void testPhiFuncs(uint32_t size, std::vector<double> const& res) {
0042   // get the list of devices on the current platform
0043   auto const& devices = cms::alpakatools::devices<Platform>();
0044   if (devices.empty()) {
0045     FAIL("No devices available for the " EDM_STRINGIZE(ALPAKA_ACCELERATOR_NAMESPACE) " backend, "
0046       "the test will be skipped.");
0047   }
0048 
0049   for (auto const& device : devices) {
0050     std::cout << "...on " << alpaka::getName(device) << "\n";
0051     Queue queue(device);
0052 
0053     auto c_h = make_host_buffer<T[]>(queue, size);
0054     alpaka::memset(queue, c_h, 0.);
0055     auto c_d = make_device_buffer<T[]>(queue, size);
0056     alpaka::memset(queue, c_d, 0.);
0057 
0058     alpaka::exec<Acc1D>(queue, WorkDiv1D{1u, 1u, 1u}, phiFuncsUnitTestsKernel(), c_d.data());
0059     alpaka::memcpy(queue, c_h, c_d);
0060     alpaka::wait(queue);
0061 
0062     constexpr T eps = 1.e-5;
0063     for (size_t i = 0; i < size; ++i) {
0064       CHECK_THAT(c_h.data()[i], Catch::Matchers::WithinAbs(res[i], eps));
0065     }
0066   }
0067 }
0068 
0069 TEST_CASE("Standard checks alpaka phi functions for the relevant data types (float and double) and for all backends") {
0070   std::vector<double> res = {0.0,
0071                              M_PI / 2.,
0072                              M_PI,
0073                              -M_PI / 2.,
0074                              M_PI / 4.,
0075                              -3. * M_PI / 4.,
0076                              -M_PI / 2.,
0077                              M_PI,
0078                              -M_PI,
0079                              -M_PI / 2.,
0080                              -M_PI / 2.,
0081                              M_PI,
0082                              -M_PI,
0083                              -M_PI / 2.};  // Expected results
0084   uint32_t size = res.size();              // Number of tests
0085 
0086   SECTION("Tests for double data type") {
0087     std::cout << "Testing phi functions for double data type...\n";
0088     testPhiFuncs<double>(size, res);
0089   }
0090 
0091   SECTION("Tests for float data type") {
0092     std::cout << "Testing phi functions for float data type...\n";
0093     testPhiFuncs<float>(size, res);
0094   }
0095 }