1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
|
#pragma once
#include <LLBM/operator.h>
namespace kernel {
template <typename OPERATOR, typename DESCRIPTOR, typename T, typename S, typename... ARGS>
__global__ void call_operator(
LatticeView<DESCRIPTOR,S> lattice
, std::size_t* cells
, std::size_t cell_count
, ARGS... args
) {
const std::size_t index = blockIdx.x * blockDim.x + threadIdx.x;
if (!(index < cell_count)) {
return;
}
const std::size_t gid = cells[index];
S f_curr[DESCRIPTOR::q];
S f_next[DESCRIPTOR::q];
S* preshifted_f[DESCRIPTOR::q];
for (unsigned iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
preshifted_f[iPop] = lattice.pop(iPop, gid);
f_curr[iPop] = *preshifted_f[iPop];
}
OPERATOR::template apply<T,S>(DESCRIPTOR(), f_curr, f_next, gid, std::forward<ARGS>(args)...);
for (unsigned iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
*preshifted_f[iPop] = f_next[iPop];
}
}
template <typename OPERATOR, typename DESCRIPTOR, typename T, typename S, typename... ARGS>
__global__ void call_operator(
LatticeView<DESCRIPTOR,S> lattice
, bool* mask
, ARGS... args
) {
const std::size_t gid = blockIdx.x * blockDim.x + threadIdx.x;
if (!(gid < lattice.cuboid.volume) || !mask[gid]) {
return;
}
S f_curr[DESCRIPTOR::q];
S f_next[DESCRIPTOR::q];
S* preshifted_f[DESCRIPTOR::q];
for (unsigned iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
preshifted_f[iPop] = lattice.pop(iPop, gid);
f_curr[iPop] = *preshifted_f[iPop];
}
OPERATOR::template apply<T,S>(DESCRIPTOR(), f_curr, f_next, gid, std::forward<ARGS>(args)...);
for (unsigned iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
*preshifted_f[iPop] = f_next[iPop];
}
}
template <typename FUNCTOR, typename DESCRIPTOR, typename T, typename S, typename... ARGS>
__global__ void call_functor(
LatticeView<DESCRIPTOR,S> lattice
, bool* mask
, ARGS... args
) {
const std::size_t gid = blockIdx.x * blockDim.x + threadIdx.x;
if (!(gid < lattice.cuboid.volume) || !mask[gid]) {
return;
}
S f_curr[DESCRIPTOR::q];
S* preshifted_f[DESCRIPTOR::q];
for (unsigned iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
preshifted_f[iPop] = lattice.pop(iPop, gid);
f_curr[iPop] = *preshifted_f[iPop];
}
FUNCTOR::template apply<T,S>(DESCRIPTOR(), f_curr, gid, std::forward<ARGS>(args)...);
}
template <typename DESCRIPTOR, typename T, typename S, typename... OPERATOR>
__global__ void call_operators(
LatticeView<DESCRIPTOR,S> lattice
, OPERATOR... ops
) {
const std::size_t gid = blockIdx.x * blockDim.x + threadIdx.x;
if (!(gid < lattice.cuboid.volume)) {
return;
}
S f_curr[DESCRIPTOR::q];
S f_next[DESCRIPTOR::q];
S* preshifted_f[DESCRIPTOR::q];
for (unsigned iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
preshifted_f[iPop] = lattice.pop(iPop, gid);
f_curr[iPop] = *preshifted_f[iPop];
}
(ops.template apply<DESCRIPTOR,T,S>(DESCRIPTOR(), f_curr, f_next, gid) || ... || false);
for (unsigned iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
*preshifted_f[iPop] = f_next[iPop];
}
}
template <typename OPERATOR, typename DESCRIPTOR, typename T, typename S, typename... ARGS>
__global__ void call_operator_using_list(
LatticeView<DESCRIPTOR,S> lattice
, std::size_t count
, ARGS... args
) {
const std::size_t index = blockIdx.x * blockDim.x + threadIdx.x;
if (!(index < count)) {
return;
}
OPERATOR::template apply<T,S>(lattice, index, count, std::forward<ARGS>(args)...);
}
template <typename OPERATOR, typename DESCRIPTOR, typename T, typename S, typename... ARGS>
__global__ void call_operator_using_list(
DESCRIPTOR descriptor
, std::size_t count
, ARGS... args
) {
const std::size_t index = blockIdx.x * blockDim.x + threadIdx.x;
if (!(index < count)) {
return;
}
OPERATOR::template apply<T,S>(descriptor, index, count, std::forward<ARGS>(args)...);
}
template <typename FUNCTOR, typename DESCRIPTOR, typename T, typename S, typename... ARGS>
__global__ void call_spatial_functor(
LatticeView<DESCRIPTOR,S> lattice
, bool* mask
, ARGS... args
) {
const std::size_t iX = blockIdx.x * blockDim.x + threadIdx.x;
const std::size_t iY = blockIdx.y * blockDim.y + threadIdx.y;
const std::size_t iZ = blockIdx.z * blockDim.z + threadIdx.z;
if (!(iX < lattice.cuboid.nX && iY < lattice.cuboid.nY && iZ < lattice.cuboid.nZ)) {
return;
}
const std::size_t gid = descriptor::gid(lattice.cuboid,iX,iY,iZ);
if (!mask[gid]) {
return;
}
S f_curr[DESCRIPTOR::q];
S* preshifted_f[DESCRIPTOR::q];
for (unsigned iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
preshifted_f[iPop] = lattice.pop(iPop, gid);
f_curr[iPop] = *preshifted_f[iPop];
}
FUNCTOR::template apply<T,S>(DESCRIPTOR(), f_curr, lattice.cuboid, gid, iX, iY, iZ, std::forward<ARGS>(args)...);
}
template <typename OPERATOR, typename DESCRIPTOR, typename T, typename S, typename... ARGS>
__global__ void call_spatial_operator(
descriptor::Cuboid<DESCRIPTOR> cuboid
, bool* mask
, ARGS... args
) {
const std::size_t iX = blockIdx.x * blockDim.x + threadIdx.x;
const std::size_t iY = blockIdx.y * blockDim.y + threadIdx.y;
const std::size_t iZ = blockIdx.z * blockDim.z + threadIdx.z;
if (!(iX < cuboid.nX && iY < cuboid.nY && iZ < cuboid.nZ)) {
return;
}
const std::size_t gid = descriptor::gid(cuboid,iX,iY,iZ);
if (!mask[gid]) {
return;
}
OPERATOR::template apply<T,S>(DESCRIPTOR(), gid, iX, iY, iZ, std::forward<ARGS>(args)...);
}
}
|