-
Notifications
You must be signed in to change notification settings - Fork 4.4k
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Let the test fail gracefully if no devies are found
If more than one device is available, run the tests on all available devices. Replace elements_with_stide with uniform_elements.
- Loading branch information
Showing
3 changed files
with
112 additions
and
103 deletions.
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
107 changes: 52 additions & 55 deletions
107
DataFormats/VertexSoA/test/alpaka/ZVertexSoA_test.dev.cc
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -1,62 +1,59 @@ | ||
#include "DataFormats/VertexSoA/interface/alpaka/ZVertexSoACollection.h" | ||
#include <alpaka/alpaka.hpp> | ||
|
||
#include "DataFormats/VertexSoA/interface/ZVertexDevice.h" | ||
#include "DataFormats/VertexSoA/interface/ZVertexHost.h" | ||
#include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h" // Check if this is really needed; code doesn't compile without it | ||
|
||
namespace ALPAKA_ACCELERATOR_NAMESPACE { | ||
using namespace alpaka; | ||
using namespace cms::alpakatools; | ||
|
||
namespace testZVertexSoAT { | ||
|
||
class TestFillKernel { | ||
public: | ||
template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>> | ||
ALPAKA_FN_ACC void operator()(TAcc const& acc, reco::ZVertexSoAView zvertex_view) const { | ||
if (cms::alpakatools::once_per_grid(acc)) { | ||
zvertex_view.nvFinal() = 420; | ||
} | ||
|
||
for (int32_t j : elements_with_stride(acc, zvertex_view.metadata().size())) { | ||
zvertex_view[j].idv() = (int16_t)j; | ||
zvertex_view[j].zv() = (float)j; | ||
zvertex_view[j].wv() = (float)j; | ||
zvertex_view[j].chi2() = (float)j; | ||
zvertex_view[j].ptv2() = (float)j; | ||
zvertex_view[j].ndof() = (int32_t)j; | ||
zvertex_view[j].sortInd() = (uint16_t)j; | ||
} | ||
#include "DataFormats/VertexSoA/interface/alpaka/ZVertexSoACollection.h" | ||
#include "HeterogeneousCore/AlpakaInterface/interface/config.h" | ||
#include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h" | ||
|
||
namespace ALPAKA_ACCELERATOR_NAMESPACE::testZVertexSoAT { | ||
|
||
class TestFillKernel { | ||
public: | ||
template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>> | ||
ALPAKA_FN_ACC void operator()(TAcc const& acc, reco::ZVertexSoAView zvertex_view) const { | ||
if (cms::alpakatools::once_per_grid(acc)) { | ||
zvertex_view.nvFinal() = 420; | ||
} | ||
}; | ||
|
||
class TestVerifyKernel { | ||
public: | ||
template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>> | ||
ALPAKA_FN_ACC void operator()(TAcc const& acc, reco::ZVertexSoAView zvertex_view) const { | ||
if (cms::alpakatools::once_per_grid(acc)) { | ||
ALPAKA_ASSERT_ACC(zvertex_view.nvFinal() == 420); | ||
} | ||
|
||
for (int32_t j : elements_with_stride(acc, zvertex_view.nvFinal())) { | ||
assert(zvertex_view[j].idv() == j); | ||
assert(zvertex_view[j].zv() - (float)j < 0.0001); | ||
assert(zvertex_view[j].wv() - (float)j < 0.0001); | ||
assert(zvertex_view[j].chi2() - (float)j < 0.0001); | ||
assert(zvertex_view[j].ptv2() - (float)j < 0.0001); | ||
assert(zvertex_view[j].ndof() == j); | ||
assert(zvertex_view[j].sortInd() == uint32_t(j)); | ||
} | ||
|
||
for (int32_t j : cms::alpakatools::uniform_elements(acc, zvertex_view.metadata().size())) { | ||
zvertex_view[j].idv() = (int16_t)j; | ||
zvertex_view[j].zv() = (float)j; | ||
zvertex_view[j].wv() = (float)j; | ||
zvertex_view[j].chi2() = (float)j; | ||
zvertex_view[j].ptv2() = (float)j; | ||
zvertex_view[j].ndof() = (int32_t)j; | ||
zvertex_view[j].sortInd() = (uint16_t)j; | ||
} | ||
} | ||
}; | ||
|
||
class TestVerifyKernel { | ||
public: | ||
template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>> | ||
ALPAKA_FN_ACC void operator()(TAcc const& acc, reco::ZVertexSoAView zvertex_view) const { | ||
if (cms::alpakatools::once_per_grid(acc)) { | ||
ALPAKA_ASSERT_ACC(zvertex_view.nvFinal() == 420); | ||
} | ||
|
||
for (int32_t j : cms::alpakatools::uniform_elements(acc, zvertex_view.nvFinal())) { | ||
assert(zvertex_view[j].idv() == j); | ||
assert(zvertex_view[j].zv() - (float)j < 0.0001); | ||
assert(zvertex_view[j].wv() - (float)j < 0.0001); | ||
assert(zvertex_view[j].chi2() - (float)j < 0.0001); | ||
assert(zvertex_view[j].ptv2() - (float)j < 0.0001); | ||
assert(zvertex_view[j].ndof() == j); | ||
assert(zvertex_view[j].sortInd() == uint32_t(j)); | ||
} | ||
}; | ||
|
||
void runKernels(reco::ZVertexSoAView zvertex_view, Queue& queue) { | ||
uint32_t items = 64; | ||
uint32_t groups = divide_up_by(zvertex_view.metadata().size(), items); | ||
auto workDiv = make_workdiv<Acc1D>(groups, items); | ||
alpaka::exec<Acc1D>(queue, workDiv, TestFillKernel{}, zvertex_view); | ||
alpaka::exec<Acc1D>(queue, workDiv, TestVerifyKernel{}, zvertex_view); | ||
} | ||
}; | ||
|
||
} // namespace testZVertexSoAT | ||
void runKernels(reco::ZVertexSoAView zvertex_view, Queue& queue) { | ||
uint32_t items = 64; | ||
uint32_t groups = cms::alpakatools::divide_up_by(zvertex_view.metadata().size(), items); | ||
auto workDiv = cms::alpakatools::make_workdiv<Acc1D>(groups, items); | ||
alpaka::exec<Acc1D>(queue, workDiv, TestFillKernel{}, zvertex_view); | ||
alpaka::exec<Acc1D>(queue, workDiv, TestVerifyKernel{}, zvertex_view); | ||
} | ||
|
||
} // namespace ALPAKA_ACCELERATOR_NAMESPACE | ||
} // namespace ALPAKA_ACCELERATOR_NAMESPACE::testZVertexSoAT |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,13 @@ | ||
#ifndef DataFormats_VertexSoA_test_alpaka_ZVertexSoA_test_h | ||
#define DataFormats_VertexSoA_test_alpaka_ZVertexSoA_test_h | ||
|
||
#include "DataFormats/VertexSoA/interface/ZVertexSoA.h" | ||
#include "HeterogeneousCore/AlpakaInterface/interface/config.h" | ||
|
||
namespace ALPAKA_ACCELERATOR_NAMESPACE::testZVertexSoAT { | ||
|
||
void runKernels(reco::ZVertexSoAView zvertex_view, Queue& queue); | ||
|
||
} // namespace ALPAKA_ACCELERATOR_NAMESPACE::testZVertexSoAT | ||
|
||
#endif // DataFormats_VertexSoA_test_alpaka_ZVertexSoA_test_h |