MaCh3  2.5.0
Reference Guide
gpuSplineUtils.cuh
Go to the documentation of this file.
1 
12 //MaCh3 included
13 #include "Manager/gpuUtils.cuh"
14 #include "Splines/SplineCommon.h"
15 
17 __host__ void SynchroniseSplines();
18 
29 __global__ void EvalOnGPU_Splines(
30  const unsigned int gpu_n_splines,
31  const short int gpu_spline_size,
32  const short int* __restrict__ gpu_paramNo_arr,
33  const unsigned int* __restrict__ gpu_nKnots_arr,
34  const float* __restrict__ gpu_coeff_many,
35  const float* __restrict__ gpu_par_val,
36  const short int* __restrict__ gpu_spline_segment,
37  float* __restrict__ gpu_weights,
38  const cudaTextureObject_t __restrict__ text_coeff_x);
39 
47 __global__ void EvalOnGPU_TF1(
48  const unsigned int gpu_n_TF1,
49  const float* __restrict__ gpu_coeffs_tf1,
50  const short int* __restrict__ gpu_paramNo_arr_tf1,
51  const float* __restrict__ gpu_par_val,
52  float* __restrict__ gpu_weights_tf1);
53 
61 __global__ void EvalOnGPU_TotWeight(
62  const int gpu_n_events,
63  const float* __restrict__ gpu_weights,
64  const float* __restrict__ gpu_weights_tf1,
65  M3::float_t* __restrict__ gpu_total_weights,
66  const cudaTextureObject_t __restrict__ text_nParamPerEvent,
67  const cudaTextureObject_t __restrict__ text_nParamPerEvent_TF1);
68 
75 {
76  public:
81  virtual ~SplineMonolithGPU();
82 
89  __host__ void InitGPU_SplineMonolith(
90  M3::float_t **cpu_total_weights,
91  int n_events,
92  unsigned int total_nknots,
93  unsigned int n_splines,
94  unsigned int n_tf1,
95  int Eve_size);
96 
114  __host__ void CopyToGPU_SplineMonolith(
115  const SplineMonoStruct* cpu_spline_handler,
116 
117  // TFI related now
118  const std::vector<float>& cpu_many_array_TF1,
119  const std::vector<short int>& cpu_paramNo_arr_TF1,
120  const int n_events,
121  const std::vector<unsigned int>& cpu_nParamPerEvent,
122  // TFI related now
123  const std::vector<unsigned int>& cpu_nParamPerEvent_TF1,
124  const int n_params,
125  const unsigned int n_splines,
126  const short int spline_size,
127  const unsigned int total_nknots,
128  const unsigned int n_tf1);
129 
132  __host__ void InitGPU_Segments(short int **segment);
133 
136  __host__ void InitGPU_Vals(float **vals);
137 
141 
157  __host__ void RunGPU_SplineMonolith(
158  M3::float_t* cpu_total_weights,
159  // Holds the changes in parameters
160  float *vals,
161  // Holds the segments for parameters
162  short int *segment);
163 
170  __host__ void CleanupPinnedMemory(M3::float_t *cpu_total_weights,
171  short int *segment, float *vals);
172 
173  private:
175  unsigned int *gpu_nParamPerEvent;
177  unsigned int *gpu_nParamPerEvent_TF1;
178 
180  float *gpu_coeff_x;
181 
184 
186  unsigned int *gpu_nKnots_arr;
187 
189  short int *gpu_paramNo_arr;
190 
194  short int *gpu_nPoints_arr;
197 
199  float* gpu_par_val;
201  short int* gpu_spline_segment;
202 
206  float *gpu_weights;
209 
211  short int cpu_spline_size;
213  unsigned int cpu_n_splines;
215  unsigned int cpu_n_TF1;
220 
221  // ******************************************
222  // TEXTURES
223  // ******************************************
225  cudaTextureObject_t text_coeff_x = 0;
227  cudaTextureObject_t text_nParamPerEvent = 0;
229  cudaTextureObject_t text_nParamPerEvent_TF1 = 0;
230 };
Contains definitions for spline coefficients and structure used in both CPU and GPU code.
Class responsible for calculating spline weight on GPU.
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...
Common CUDA utilities and definitions for shared GPU functionality.
double float_t
Definition: Core.h:37
KS: Struct storing information for spline monolith.
Definition: SplineCommon.h:35