/
Backend_par.hpp
212 lines (180 loc) · 5.54 KB
/
Backend_par.hpp
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
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
#pragma once
#ifndef SPIRIT_CORE_ENGINE_BACKEND_PAR_HPP
#define SPIRIT_CORE_ENGINE_BACKEND_PAR_HPP
#include <engine/Vectormath_Defines.hpp>
// clang-format off
#ifdef SPIRIT_USE_CUDA
#define THRUST_IGNORE_CUB_VERSION_CHECK
#include <cub/cub.cuh>
#define SPIRIT_LAMBDA __device__
#else
#define SPIRIT_LAMBDA
#endif
// clang-format on
namespace Spirit::Engine::Backend::par
{
#ifdef SPIRIT_USE_CUDA
template<typename F>
__global__ void cu_apply( int N, F f )
{
int idx = blockIdx.x * blockDim.x + threadIdx.x;
if( idx < N )
f( idx );
}
// vf1[i] = f( vf2[i] )
template<typename F>
void apply( int N, F f )
{
cu_apply<<<( N + 1023 ) / 1024, 1024>>>( N, f );
CU_CHECK_AND_SYNC();
}
template<typename F>
scalar reduce( int N, const F f )
{
static scalarfield sf( N, 0 );
// Vectormath::fill(sf, 0);
if( sf.size() != N )
sf.resize( N );
auto s = sf.data();
apply( N, [f, s] SPIRIT_LAMBDA( int idx ) { s[idx] = f( idx ); } );
static scalarfield ret( 1, 0 );
// Determine temporary storage size and allocate
void * d_temp_storage = NULL;
size_t temp_storage_bytes = 0;
cub::DeviceReduce::Sum( d_temp_storage, temp_storage_bytes, sf.data(), ret.data(), sf.size() );
cudaMalloc( &d_temp_storage, temp_storage_bytes );
// Reduction
cub::DeviceReduce::Sum( d_temp_storage, temp_storage_bytes, sf.data(), ret.data(), sf.size() );
CU_CHECK_AND_SYNC();
cudaFree( d_temp_storage );
return ret[0];
}
template<typename A, typename F>
scalar reduce( const field<A> & vf1, const F f )
{
// TODO: remove the reliance on a temporary scalar field (maybe thrust::dot with generalized operations)
// We also use this workaround for a single field as argument, because cub does not support non-commutative
// reduction operations
int n = vf1.size();
static scalarfield sf( n, 0 );
// Vectormath::fill(sf, 0);
if( sf.size() != vf1.size() )
sf.resize( vf1.size() );
auto s = sf.data();
auto v1 = vf1.data();
apply( n, [f, s, v1] SPIRIT_LAMBDA( int idx ) { s[idx] = f( v1[idx] ); } );
static scalarfield ret( 1, 0 );
// Vectormath::fill(ret, 0);
// Determine temporary storage size and allocate
void * d_temp_storage = NULL;
size_t temp_storage_bytes = 0;
cub::DeviceReduce::Sum( d_temp_storage, temp_storage_bytes, sf.data(), ret.data(), sf.size() );
cudaMalloc( &d_temp_storage, temp_storage_bytes );
// Reduction
cub::DeviceReduce::Sum( d_temp_storage, temp_storage_bytes, sf.data(), ret.data(), sf.size() );
CU_CHECK_AND_SYNC();
cudaFree( d_temp_storage );
return ret[0];
}
template<typename A, typename B, typename F>
scalar reduce( const field<A> & vf1, const field<B> & vf2, const F f )
{
// TODO: remove the reliance on a temporary scalar field (maybe thrust::dot with generalized operations)
int n = vf1.size();
static scalarfield sf( n, 0 );
// Vectormath::fill(sf, 0);
if( sf.size() != vf1.size() )
sf.resize( vf1.size() );
auto s = sf.data();
auto v1 = vf1.data();
auto v2 = vf2.data();
apply( n, [f, s, v1, v2] SPIRIT_LAMBDA( int idx ) { s[idx] = f( v1[idx], v2[idx] ); } );
static scalarfield ret( 1, 0 );
// Vectormath::fill(ret, 0);
// Determine temporary storage size and allocate
void * d_temp_storage = NULL;
size_t temp_storage_bytes = 0;
cub::DeviceReduce::Sum( d_temp_storage, temp_storage_bytes, sf.data(), ret.data(), sf.size() );
cudaMalloc( &d_temp_storage, temp_storage_bytes );
// Reduction
cub::DeviceReduce::Sum( d_temp_storage, temp_storage_bytes, sf.data(), ret.data(), sf.size() );
CU_CHECK_AND_SYNC();
cudaFree( d_temp_storage );
return ret[0];
}
// vf1[i] = f( vf2[i] )
template<typename A, typename B, typename F>
__global__ void cu_set_lambda( A * vf1, const B * vf2, F f, int N )
{
int idx = blockIdx.x * blockDim.x + threadIdx.x;
if( idx < N )
{
vf1[idx] = f( vf2[idx] );
}
}
// vf1[i] = f( vf2[i] )
template<typename A, typename B, typename F>
void set( field<A> & vf1, const field<B> & vf2, F f )
{
int N = vf1.size();
cu_set_lambda<<<( N + 1023 ) / 1024, 1024>>>( vf1.data(), vf2.data(), f, N );
CU_CHECK_AND_SYNC();
}
#else
template<typename F>
scalar reduce( int N, const F f )
{
scalar res = 0;
#pragma omp parallel for reduction( + : res )
for( unsigned int idx = 0; idx < N; ++idx )
{
res += f( idx );
}
return res;
}
template<typename A, typename F>
scalar reduce( const field<A> & vf1, const F f )
{
scalar res = 0;
#pragma omp parallel for reduction( + : res )
for( unsigned int idx = 0; idx < vf1.size(); ++idx )
{
res += f( vf1[idx] );
}
return res;
}
// result = sum_i f( vf1[i], vf2[i] )
template<typename A, typename B, typename F>
scalar reduce( const field<A> & vf1, const field<B> & vf2, const F & f )
{
scalar res = 0;
#pragma omp parallel for reduction( + : res )
for( unsigned int idx = 0; idx < vf1.size(); ++idx )
{
res += f( vf1[idx], vf2[idx] );
}
return res;
}
// vf1[i] = f( vf2[i] )
template<typename A, typename B, typename F>
void set( field<A> & vf1, const field<B> & vf2, const F & f )
{
#pragma omp parallel for
for( unsigned int idx = 0; idx < vf1.size(); ++idx )
{
vf1[idx] = f( vf2[idx] );
}
}
// f( vf1[idx], idx ) for all i
template<typename F>
void apply( int N, const F & f )
{
#pragma omp parallel for
for( unsigned int idx = 0; idx < N; ++idx )
{
f( idx );
}
}
#endif
} // namespace Engine::Backend::par
#endif