src
HeterogeneousCore
AlpakaTest
plugins
alpaka
TestAlgo.dev.cc
Go to the documentation of this file.
1
// Check that ALPAKA_HOST_ONLY is not defined during device compilation:
2
#ifdef ALPAKA_HOST_ONLY
3
#error ALPAKA_HOST_ONLY defined in device compilation
4
#endif
5
6
#include <alpaka/alpaka.hpp>
7
8
#include "
DataFormats/PortableTestObjects/interface/alpaka/TestDeviceCollection.h
"
9
#include "
HeterogeneousCore/AlpakaInterface/interface/config.h
"
10
#include "
HeterogeneousCore/AlpakaInterface/interface/traits.h
"
11
#include "
HeterogeneousCore/AlpakaInterface/interface/workdivision.h
"
12
13
#include "
TestAlgo.h
"
14
15
namespace
ALPAKA_ACCELERATOR_NAMESPACE
{
16
17
using namespace
cms::alpakatools
;
18
19
class
TestAlgoKernel
{
20
public
:
21
template
<
typename
TAcc,
typename
= std::enable_if_t<alpaka::isAccelerator<TAcc>>>
22
ALPAKA_FN_ACC
void
operator()
(TAcc
const
& acc,
23
portabletest::TestDeviceCollection::View
view
,
24
int32_t
size
,
25
double
xvalue)
const
{
26
// global index of the thread within the grid
27
const
int32_t thread = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc)[0u];
28
const
portabletest::Matrix
matrix
{{1, 2, 3, 4, 5, 6}, {2, 4, 6, 8, 10, 12}, {3, 6, 9, 12, 15, 18}};
29
30
// set this only once in the whole kernel grid
31
if
(thread == 0) {
32
view
.r() = 1.;
33
}
34
35
// make a strided loop over the kernel grid, covering up to "size" elements
36
for
(int32_t
i
:
elements_with_stride
(acc,
size
)) {
37
view
[
i
] = {xvalue, 0., 0.,
i
,
matrix
*
i
};
38
}
39
}
40
};
41
42
void
TestAlgo::fill
(Queue&
queue
,
portabletest::TestDeviceCollection
&
collection
,
double
xvalue)
const
{
43
// use 64 items per group (this value is arbitrary, but it's a reasonable starting point)
44
uint32_t
items
= 64;
45
46
// use as many groups as needed to cover the whole problem
47
uint32_t groups =
divide_up_by
(
collection
->metadata().size(),
items
);
48
49
// map items to
50
// - threads with a single element per thread on a GPU backend
51
// - elements within a single thread on a CPU backend
52
auto
workDiv = make_workdiv<Acc1D>(groups,
items
);
53
54
alpaka::exec<Acc1D>(
queue
, workDiv,
TestAlgoKernel
{},
collection
.view(),
collection
->metadata().size(), xvalue);
55
}
56
57
}
// namespace ALPAKA_ACCELERATOR_NAMESPACE
findQualityFiles.size
size
Write out results.
Definition:
findQualityFiles.py:443
TestDeviceCollection.h
makeMuonMisalignmentScenario.matrix
list matrix
Definition:
makeMuonMisalignmentScenario.py:141
mps_fire.i
i
Definition:
mps_fire.py:429
cms::alpakatools::divide_up_by
constexpr Idx divide_up_by(Idx value, Idx divisor)
Definition:
workdivision.h:19
portabletest::Matrix
Eigen::Matrix< double, 3, 6 > Matrix
Definition:
TestSoA.h:13
pv::view
view
Definition:
PreparePVTrends.h:58
workdivision.h
ALPAKA_ACCELERATOR_NAMESPACE::TestAlgoKernel::operator()
ALPAKA_FN_ACC void operator()(TAcc const &acc, portabletest::TestDeviceCollection::View view, int32_t size, double xvalue) const
Definition:
TestAlgo.dev.cc:22
createBeamHaloJobs.queue
queue
Definition:
createBeamHaloJobs.py:343
cms::alpakatools
Definition:
PortableCollection.h:44
ALPAKA_ACCELERATOR_NAMESPACE
Definition:
PortableCollection.h:15
PortableDeviceCollection::View
typename Layout::View View
Definition:
PortableDeviceCollection.h:21
universalConfigTemplate.collection
collection
Definition:
universalConfigTemplate.py:81
traits.h
ALPAKA_ACCELERATOR_NAMESPACE::TestAlgo::fill
void fill(Queue &queue, portabletest::TestDeviceCollection &collection, double xvalue=0.) const
Definition:
TestAlgo.dev.cc:42
mps_monitormerge.items
list items
Definition:
mps_monitormerge.py:29
cms::alpakatools::elements_with_stride
Definition:
workdivision.h:79
TestAlgo.h
PortableDeviceCollection
Definition:
PortableDeviceCollection.h:15
ALPAKA_ACCELERATOR_NAMESPACE::TestAlgoKernel
Definition:
TestAlgo.dev.cc:19
config.h
Generated for CMSSW Reference Manual by
1.8.14