5 #if defined(__CUDA_ARCH__)
6 #if __CUDA_ARCH__ >= 1200
7 #pragma message("Compiling with CUDA Architecture: 12.x")
8 #elif __CUDA_ARCH__ >= 1100
9 #pragma message("Compiling with CUDA Architecture: 11.x")
10 #elif __CUDA_ARCH__ >= 1000
11 #pragma message("Compiling with CUDA Architecture: 10.x")
12 #elif __CUDA_ARCH__ >= 900
13 #pragma message("Compiling with CUDA Architecture: 9.x")
14 #elif __CUDA_ARCH__ >= 800
15 #pragma message("Compiling with CUDA Architecture: 8.x")
16 #elif __CUDA_ARCH__ >= 750
17 #pragma message("Compiling with CUDA Architecture: 7.5")
18 #elif __CUDA_ARCH__ >= 730
19 #pragma message("Compiling with CUDA Architecture: 7.3")
20 #elif __CUDA_ARCH__ >= 720
21 #pragma message("Compiling with CUDA Architecture: 7.2")
22 #elif __CUDA_ARCH__ >= 710
23 #pragma message("Compiling with CUDA Architecture: 7.1")
24 #elif __CUDA_ARCH__ >= 700
25 #pragma message("Compiling with CUDA Architecture: 7.x")
26 #elif __CUDA_ARCH__ >= 650
27 #pragma message("Compiling with CUDA Architecture: 6.5")
28 #elif __CUDA_ARCH__ >= 600
29 #pragma message("Compiling with CUDA Architecture: 6.x")
30 #elif __CUDA_ARCH__ >= 530
31 #pragma message("Compiling with CUDA Architecture: 5.3")
32 #elif __CUDA_ARCH__ >= 520
33 #pragma message("Compiling with CUDA Architecture: 5.2")
34 #elif __CUDA_ARCH__ >= 510
35 #pragma message("Compiling with CUDA Architecture: 5.1")
36 #elif __CUDA_ARCH__ >= 500
37 #pragma message("Compiling with CUDA Architecture: 5.x")
38 #elif __CUDA_ARCH__ >= 400
39 #pragma message("Compiling with CUDA Architecture: 4.x")
40 #elif __CUDA_ARCH__ >= 300
41 #pragma message("Compiling with CUDA Architecture: 3.x")
43 #pragma message("Compiling with CUDA Architecture: < 3.x")
51 cudaDeviceSynchronize();
113 unsigned int total_nknots,
114 unsigned int n_splines,
122 cudaMalloc((
void **) &
gpu_nKnots_arr, n_splines*
sizeof(
unsigned int));
125 cudaMalloc((
void **) &
gpu_coeff_x, Eve_size*
sizeof(
float));
132 cudaMalloc((
void **) &
gpu_weights, n_splines*
sizeof(
float));
146 cudaMallocHost((
void **) cpu_total_weights, n_events*
sizeof(
M3::float_t));
162 printf(
"Allocated %i entries for paramNo and nKnots arrays, size = %f MB\n",
163 n_splines,
static_cast<double>(
sizeof(
short int) * n_splines +
sizeof(
unsigned int) * n_splines) / 1.0e6);
164 printf(
"Allocated %i entries for x coeff arrays, size = %f MB\n",
165 Eve_size,
static_cast<double>(
sizeof(
float) * Eve_size) / 1.0e6);
166 printf(
"Allocated %i entries for {ybcd} coeff arrays, size = %f MB\n",
167 _nCoeff_ * total_nknots,
static_cast<double>(
sizeof(
float) *
_nCoeff_ * total_nknots) / 1.0e6);
168 printf(
"Allocated %i entries for TF1 coefficient arrays, size = %f MB\n",
181 cudaMallocHost((
void **) segment,
cpu_n_params*
sizeof(
short int));
190 cudaMallocHost((
void **) vals,
cpu_n_params*
sizeof(
float));
204 const std::vector<float>& cpu_many_array_TF1,
205 const std::vector<short int>& cpu_paramNo_arr_TF1,
207 const std::vector<unsigned int>& cpu_nParamPerEvent,
209 const std::vector<unsigned int>& cpu_nParamPerEvent_TF1,
212 const unsigned int n_splines,
213 const short int spline_size,
214 const unsigned int total_nknots,
215 const unsigned int n_tf1) {
229 cudaMalloc(&
gpu_par_val, n_params *
sizeof(
float));
236 cudaMemcpy(
gpu_coeff_x, cpu_spline_handler->
coeff_x.data(),
sizeof(
float)*spline_size*n_params, cudaMemcpyHostToDevice);
241 struct cudaResourceDesc resDesc_coeff_x;
242 memset(&resDesc_coeff_x, 0,
sizeof(resDesc_coeff_x));
243 resDesc_coeff_x.resType = cudaResourceTypeLinear;
245 resDesc_coeff_x.res.linear.desc = cudaCreateChannelDesc<float>();
246 resDesc_coeff_x.res.linear.sizeInBytes =
sizeof(float)*spline_size*n_params;
249 struct cudaTextureDesc texDesc_coeff_x;
250 memset(&texDesc_coeff_x, 0,
sizeof(texDesc_coeff_x));
251 texDesc_coeff_x.readMode = cudaReadModeElementType;
254 cudaCreateTextureObject(&
text_coeff_x, &resDesc_coeff_x, &texDesc_coeff_x,
nullptr);
262 cudaMemcpy(
gpu_nKnots_arr, cpu_spline_handler->
nKnots_arr.data(), n_splines*
sizeof(
unsigned int), cudaMemcpyHostToDevice);
271 cudaMemcpy(
gpu_paramNo_TF1_arr, cpu_paramNo_arr_TF1.data(), n_tf1*
sizeof(
short int), cudaMemcpyHostToDevice);
275 cudaMemcpy(
gpu_nParamPerEvent, cpu_nParamPerEvent.data(), 2*n_events*
sizeof(
unsigned int), cudaMemcpyHostToDevice);
280 struct cudaResourceDesc resDesc_nParamPerEvent;
281 memset(&resDesc_nParamPerEvent, 0,
sizeof(resDesc_nParamPerEvent));
282 resDesc_nParamPerEvent.resType = cudaResourceTypeLinear;
284 resDesc_nParamPerEvent.res.linear.desc = cudaCreateChannelDesc<unsigned int>();
285 resDesc_nParamPerEvent.res.linear.sizeInBytes = 2*n_events*
sizeof(
unsigned int);
288 struct cudaTextureDesc texDesc_nParamPerEvent;
289 memset(&texDesc_nParamPerEvent, 0,
sizeof(texDesc_nParamPerEvent));
290 texDesc_nParamPerEvent.readMode = cudaReadModeElementType;
293 cudaCreateTextureObject(&
text_nParamPerEvent, &resDesc_nParamPerEvent, &texDesc_nParamPerEvent,
nullptr);
297 cudaMemcpy(
gpu_nParamPerEvent_TF1, cpu_nParamPerEvent_TF1.data(), 2*n_events*
sizeof(
unsigned int), cudaMemcpyHostToDevice);
302 struct cudaResourceDesc resDesc_nParamPerEvent_tf1;
303 memset(&resDesc_nParamPerEvent_tf1, 0,
sizeof(resDesc_nParamPerEvent_tf1));
304 resDesc_nParamPerEvent_tf1.resType = cudaResourceTypeLinear;
306 resDesc_nParamPerEvent_tf1.res.linear.desc = cudaCreateChannelDesc<unsigned int>();
307 resDesc_nParamPerEvent_tf1.res.linear.sizeInBytes = 2*n_events*
sizeof(
unsigned int);
310 struct cudaTextureDesc texDesc_nParamPerEvent_tf1;
311 memset(&texDesc_nParamPerEvent_tf1, 0,
sizeof(texDesc_nParamPerEvent_tf1));
312 texDesc_nParamPerEvent_tf1.readMode = cudaReadModeElementType;
315 cudaCreateTextureObject(&
text_nParamPerEvent_TF1, &resDesc_nParamPerEvent_tf1, &texDesc_nParamPerEvent_tf1,
nullptr);
331 const unsigned int gpu_n_splines,
332 const short int gpu_spline_size,
333 const short int* __restrict__ gpu_paramNo_arr,
334 const unsigned int* __restrict__ gpu_nKnots_arr,
335 const float* __restrict__ gpu_coeff_many,
336 const float* __restrict__ gpu_par_val,
337 const short int* __restrict__ gpu_spline_segment,
338 float* __restrict__ gpu_weights,
339 const cudaTextureObject_t __restrict__ text_coeff_x) {
342 const unsigned int splineNum = (blockIdx.x * blockDim.x + threadIdx.x);
345 if (splineNum < gpu_n_splines) {
350 const short int Param = gpu_paramNo_arr[splineNum];
353 const short int segment = gpu_spline_segment[Param];
356 const short int segment_X = Param*gpu_spline_size+segment;
359 const unsigned int CurrentKnotPos = gpu_nKnots_arr[splineNum]*
_nCoeff_+segment*
_nCoeff_;
363 const float fY = gpu_coeff_many[CurrentKnotPos];
364 const float fB = gpu_coeff_many[CurrentKnotPos + 1];
365 const float fC = gpu_coeff_many[CurrentKnotPos + 2];
366 const float fD = gpu_coeff_many[CurrentKnotPos + 3];
368 const float dx = gpu_par_val[Param] - tex1Dfetch<float>(text_coeff_x, segment_X);
371 gpu_weights[splineNum] = fmaf(dx, fmaf(dx, fmaf(dx, fD, fC), fB), fY);
384 const unsigned int gpu_n_TF1,
385 const float* __restrict__ gpu_coeffs_tf1,
386 const short int* __restrict__ gpu_paramNo_arr_tf1,
387 const float* __restrict__ gpu_par_val,
388 float* __restrict__ gpu_weights_tf1) {
391 const unsigned int tf1Num = (blockIdx.x * blockDim.x + threadIdx.x);
393 if (tf1Num < gpu_n_TF1) {
395 const float x = gpu_par_val[gpu_paramNo_arr_tf1[tf1Num]];
398 const unsigned int TF1_Index = tf1Num *
_nTF1Coeff_;
399 const float a = gpu_coeffs_tf1[TF1_Index];
400 const float b = gpu_coeffs_tf1[TF1_Index+1];
402 gpu_weights_tf1[tf1Num] = fmaf(a, x, b);
411 const int gpu_n_events,
412 const float* __restrict__ gpu_weights,
413 const float* __restrict__ gpu_weights_tf1,
417 const cudaTextureObject_t __restrict__ text_nParamPerEvent,
418 const cudaTextureObject_t __restrict__ text_nParamPerEvent_TF1) {
420 const unsigned int EventNum = (blockIdx.x * blockDim.x + threadIdx.x);
422 if(EventNum < gpu_n_events)
424 float local_total_weight = 1.f;
426 const unsigned int EventOffset = 2 * EventNum;
428 for (
unsigned int id = 0; id < tex1Dfetch<unsigned int>(text_nParamPerEvent, EventOffset); ++id) {
429 local_total_weight *= gpu_weights[tex1Dfetch<unsigned int>(text_nParamPerEvent, EventOffset+1) + id];
432 for (
unsigned int id = 0; id < tex1Dfetch<unsigned int>(text_nParamPerEvent_TF1, EventOffset); ++id) {
433 local_total_weight *= gpu_weights_tf1[tex1Dfetch<unsigned int>(text_nParamPerEvent_TF1, EventOffset+1) + id];
435 gpu_total_weights[EventNum] =
static_cast<M3::float_t>(local_total_weight);
448 short int *segment) {
467 EvalOnGPU_Splines<<<grid_size, block_size>>>(
481 grid_size.x = (
cpu_n_TF1 / block_size.x) + 1;
482 EvalOnGPU_TF1<<<grid_size, block_size>>>(
493 EvalOnGPU_TotWeight<<<grid_size, block_size>>>(
511 printf(
"Copied GPU total weights to CPU with SUCCESS (drink more tea)\n");
512 printf(
"Released calculated response from GPU with SUCCESS (drink most tea)\n");
519 short int *segment,
float *vals) {
521 cudaFreeHost(cpu_total_weights);
522 cudaFreeHost(segment);
525 cpu_total_weights =
nullptr;
constexpr int _nCoeff_
KS: We store coefficients {y,b,c,d} in one array one by one, this is only to define it once rather th...
constexpr int _nTF1Coeff_
KS: For TF1 we store at most 5 coefficients, we could make it more flexible but for now define it her...
float * gpu_weights
GPU arrays to hold weight for each spline.
unsigned int * gpu_nParamPerEvent
KS: GPU map keeping track how many parameters applies to each event, we keep two numbers here {number...
__host__ void RunGPU_SplineMonolith(M3::float_t *cpu_total_weights, float *vals, short int *segment)
Run the GPU code for the separate many arrays. As in separate {x}, {y,b,c,d} arrays Pass the segment ...
unsigned int cpu_n_TF1
Number of tf1 living on CPU.
virtual ~SplineMonolithGPU()
This function deallocates the resources allocated for the separate {x} and {ybcd} arrays in the and T...
unsigned int * gpu_nParamPerEvent_TF1
KS: GPU map keeping track how many parameters applies to each event, we keep two numbers here {number...
float * gpu_weights_tf1
GPU arrays to hold weight for each TF1.
unsigned int * gpu_nKnots_arr
KS: GPU Number of knots per spline.
float * gpu_par_val
CW: parameter value on GPU.
float * gpu_coeff_TF1_many
GPU arrays to hold TF1 coefficients.
__host__ void InitGPU_Vals(float **vals)
Allocate memory for spline segments.
cudaTextureObject_t text_nParamPerEvent_TF1
KS: Map keeping track how many parameters applies to each event, we keep two numbers here {number of ...
float * gpu_coeff_many
GPU arrays to hold other coefficients.
__host__ void CleanupPinnedMemory(M3::float_t *cpu_total_weights, short int *segment, float *vals)
Clean up pinned variables at CPU.
short int * gpu_paramNo_TF1_arr
CW: GPU array with the number of points per TF1 object.
__host__ void InitGPU_SplineMonolith(M3::float_t **cpu_total_weights, int n_events, unsigned int total_nknots, unsigned int n_splines, unsigned int n_tf1, int Eve_size)
Allocate memory on gpu for spline monolith.
cudaTextureObject_t text_coeff_x
KS: Textures are L1 cache variables which are well optimised for fetching. Make texture only for vari...
int cpu_n_params
Number of params living on CPU.
__host__ void CopyToGPU_SplineMonolith(const SplineMonoStruct *cpu_spline_handler, const std::vector< float > &cpu_many_array_TF1, const std::vector< short int > &cpu_paramNo_arr_TF1, const int n_events, const std::vector< unsigned int > &cpu_nParamPerEvent, const std::vector< unsigned int > &cpu_nParamPerEvent_TF1, const int n_params, const unsigned int n_splines, const short int spline_size, const unsigned int total_nknots, const unsigned int n_tf1)
Copies data from CPU to GPU for the spline monolith.
__host__ void InitGPU_Segments(short int **segment)
Allocate memory for spline segments.
int cpu_n_events
Number of events living on CPU.
short int * gpu_nPoints_arr
GPU arrays to hold number of points.
unsigned int cpu_n_splines
Number of splines living on CPU.
SplineMonolithGPU()
constructor
short int * gpu_paramNo_arr
CW: GPU array with the number of points per spline (not per spline point!)
float * gpu_coeff_x
KS: GPU arrays to hold X coefficient.
cudaTextureObject_t text_nParamPerEvent
KS: Map keeping track how many parameters applies to each event, we keep two numbers here {number of ...
short int cpu_spline_size
Size of splines living on CPU.
M3::float_t * gpu_total_weights
GPU arrays to hold weight for event.
short int * gpu_spline_segment
CW: Spline segment on GPU.
__global__ void EvalOnGPU_TF1(const unsigned int gpu_n_TF1, const float *__restrict__ gpu_coeffs_tf1, const short int *__restrict__ gpu_paramNo_arr_tf1, const float *__restrict__ gpu_par_val, float *__restrict__ gpu_weights_tf1)
Evaluate the TF1 on the GPU Using first order polynomial Polynomial form: w(x) = a0 + a1·x.
__host__ void SynchroniseSplines()
Make sure all Cuda threads finished execution.
__global__ void EvalOnGPU_TotWeight(const int gpu_n_events, const float *__restrict__ gpu_weights, const float *__restrict__ gpu_weights_tf1, M3::float_t *__restrict__ gpu_total_weights, const cudaTextureObject_t __restrict__ text_nParamPerEvent, const cudaTextureObject_t __restrict__ text_nParamPerEvent_TF1)
KS: Evaluate the total spline event weight on the GPU, as in most cases GPU is faster,...
__global__ void EvalOnGPU_Splines(const unsigned int gpu_n_splines, const short int gpu_spline_size, const short int *__restrict__ gpu_paramNo_arr, const unsigned int *__restrict__ gpu_nKnots_arr, const float *__restrict__ gpu_coeff_many, const float *__restrict__ gpu_par_val, const short int *__restrict__ gpu_spline_segment, float *__restrict__ gpu_weights, const cudaTextureObject_t __restrict__ text_coeff_x)
Evaluate the spline on the GPU Using one {y,b,c,d} array and one {x} array Should be most efficient a...
MaCh3 event-by-event cross-section spline code.
void checkGpuMem()
KS: Get some fancy info about VRAM usage.
void PrintNdevices()
KS: Get some fancy info about GPU.
#define _BlockSize_
KS: Need it for shared memory, there is way to use dynamic shared memory but I am lazy right now.
KS: Struct storing information for spline monolith.
std::vector< unsigned int > nKnots_arr
KS: CPU Number of knots per spline.
std::vector< float > coeff_x
KS: CPU arrays to hold X coefficient.
std::vector< float > coeff_many
CPU arrays to hold other coefficients.
std::vector< short int > paramNo_arr
CW: CPU array with the number of points per spline (not per spline point!)