![]() |
MaCh3
2.5.0
Reference Guide
|
Class responsible for calculating spline weight on GPU. More...
Public Member Functions | |
| SplineMonolithGPU () | |
| constructor More... | |
| virtual | ~SplineMonolithGPU () |
| This function deallocates the resources allocated for the separate {x} and {ybcd} arrays in the and TF1 stuff at GPU. More... | |
| __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. More... | |
| __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. More... | |
| __host__ void | InitGPU_Segments (short int **segment) |
| Allocate memory for spline segments. More... | |
| __host__ void | InitGPU_Vals (float **vals) |
| Allocate memory for spline segments. More... | |
| __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 and the parameter values (binary search already performed in SplineBase::FindSplineSegment() More... | |
| __host__ void | CleanupPinnedMemory (M3::float_t *cpu_total_weights, short int *segment, float *vals) |
| Clean up pinned variables at CPU. More... | |
Private Attributes | |
| unsigned int * | gpu_nParamPerEvent |
| KS: GPU map keeping track how many parameters applies to each event, we keep two numbers here {number of splines per event, index where splines start for a given event}. More... | |
| unsigned int * | gpu_nParamPerEvent_TF1 |
| KS: GPU map keeping track how many parameters applies to each event, we keep two numbers here {number of TF1 per event, index where TF1 start for a given event}. More... | |
| float * | gpu_coeff_x |
| KS: GPU arrays to hold X coefficient. More... | |
| float * | gpu_coeff_many |
| GPU arrays to hold other coefficients. More... | |
| unsigned int * | gpu_nKnots_arr |
| KS: GPU Number of knots per spline. More... | |
| short int * | gpu_paramNo_arr |
| CW: GPU array with the number of points per spline (not per spline point!) More... | |
| float * | gpu_coeff_TF1_many |
| GPU arrays to hold TF1 coefficients. More... | |
| short int * | gpu_nPoints_arr |
| GPU arrays to hold number of points. More... | |
| short int * | gpu_paramNo_TF1_arr |
| CW: GPU array with the number of points per TF1 object. More... | |
| float * | gpu_par_val |
| CW: parameter value on GPU. More... | |
| short int * | gpu_spline_segment |
| CW: Spline segment on GPU. More... | |
| M3::float_t * | gpu_total_weights |
| GPU arrays to hold weight for event. More... | |
| float * | gpu_weights |
| GPU arrays to hold weight for each spline. More... | |
| float * | gpu_weights_tf1 |
| GPU arrays to hold weight for each TF1. More... | |
| short int | cpu_spline_size |
| Size of splines living on CPU. More... | |
| unsigned int | cpu_n_splines |
| Number of splines living on CPU. More... | |
| unsigned int | cpu_n_TF1 |
| Number of tf1 living on CPU. More... | |
| int | cpu_n_params |
| Number of params living on CPU. More... | |
| int | cpu_n_events |
| Number of events living on CPU. More... | |
| cudaTextureObject_t | text_coeff_x = 0 |
| KS: Textures are L1 cache variables which are well optimised for fetching. Make texture only for variables you often access but rarely overwrite. There are limits on texture memory so don't use huge arrays. More... | |
| cudaTextureObject_t | text_nParamPerEvent = 0 |
| KS: Map keeping track how many parameters applies to each event, we keep two numbers here {number of splines per event, index where splines start for a given event}. More... | |
| cudaTextureObject_t | text_nParamPerEvent_TF1 = 0 |
| KS: Map keeping track how many parameters applies to each event, we keep two numbers here {number of TF1 per event, index where TF1 start for a given event}. More... | |
Class responsible for calculating spline weight on GPU.
Definition at line 74 of file gpuSplineUtils.cuh.
| SplineMonolithGPU::SplineMonolithGPU | ( | ) |
constructor
Number of events living on CPU
Definition at line 58 of file gpuSplineUtils.cu.
|
virtual |
This function deallocates the resources allocated for the separate {x} and {ybcd} arrays in the and TF1 stuff at GPU.
Definition at line 81 of file gpuSplineUtils.cu.
| __host__ void SplineMonolithGPU::CleanupPinnedMemory | ( | M3::float_t * | cpu_total_weights, |
| short int * | segment, | ||
| float * | vals | ||
| ) |
Clean up pinned variables at CPU.
| cpu_total_weights | Pointer to the total weights array on the CPU (used if Weight_On_SplineBySpline_Basis is not defined). |
| segment | Found spline segment for each parameter |
| vals | Value to which we want reweight for each parameter |
Definition at line 518 of file gpuSplineUtils.cu.
| __host__ void SplineMonolithGPU::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.
This function transfers the necessary spline data and parameters from the CPU to the GPU, including TF1-related arrays and parameters. This setup is crucial for spline evaluations on the GPU.
| cpu_spline_handler | Pointer to the structure managing spline data on the CPU. |
| cpu_many_array_TF1 | Array of TF1 parameters on the CPU. |
| cpu_paramNo_arr_TF1 | Array containing parameter numbers for TF1 objects. |
| n_events | Number of events, necessary for correct data handling. |
| cpu_nParamPerEvent | Array indicating the number of parameters per event. |
| cpu_nParamPerEvent_TF1 | Array indicating the number of parameters per TF1 object. |
| n_params | Total number of parameters across all splines. |
| n_splines | Total number of spline objects. |
| spline_size | Size of each spline object. |
| total_nknots | Total number of knots across all splines. |
| n_tf1 | Total number of TF1 objects. |
Size of splines living
Definition at line 200 of file gpuSplineUtils.cu.
| __host__ void SplineMonolithGPU::InitGPU_Segments | ( | short int ** | segment | ) |
Allocate memory for spline segments.
| segment | Found spline segment for each parameter |
Definition at line 178 of file gpuSplineUtils.cu.
| __host__ void SplineMonolithGPU::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.
| cpu_total_weights | KS: Rather than allocate memory in standard way this fancy cuda tool allows to pin host memory which make memory transfer faster |
| n_events | Number of events, this is necessary to allocate correctly memory |
| total_nknots | Total number of knots in all splines summed |
| n_splines | Total number of spline objects, not knots |
| n_tf1 | Total number of TF1 objects, not coefficients |
Definition at line 110 of file gpuSplineUtils.cu.
| __host__ void SplineMonolithGPU::InitGPU_Vals | ( | float ** | vals | ) |
Allocate memory for spline segments.
| vals | Value to which we want reweight for each parameter |
Definition at line 187 of file gpuSplineUtils.cu.
| __host__ void SplineMonolithGPU::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 and the parameter values (binary search already performed in SplineBase::FindSplineSegment()
Executes the GPU code for calculating spline weights.
This function runs the GPU computation for the spline monolith. It assumes that the appropriate segment has already been identified through binary search in the SplineBase::FindSplineSegment() function.
| cpu_weights | Pointer to the array of weights on the CPU (used if Weight_On_SplineBySpline_Basis is defined). |
| cpu_weights_tf1 | Pointer to the array of TF1 weights (used if Weight_On_SplineBySpline_Basis is defined). |
| cpu_total_weights | Pointer to the total weights array (used if Weight_On_SplineBySpline_Basis is not defined). |
| vals | Pointer to an array holding the parameter values to be processed. |
| segment | Pointer to an array containing segment indices for parameters. |
Definition at line 443 of file gpuSplineUtils.cu.
|
private |
Number of events living on CPU.
Definition at line 219 of file gpuSplineUtils.cuh.
|
private |
Number of params living on CPU.
Definition at line 217 of file gpuSplineUtils.cuh.
|
private |
Number of splines living on CPU.
Definition at line 213 of file gpuSplineUtils.cuh.
|
private |
Number of tf1 living on CPU.
Definition at line 215 of file gpuSplineUtils.cuh.
|
private |
Size of splines living on CPU.
Definition at line 211 of file gpuSplineUtils.cuh.
|
private |
GPU arrays to hold other coefficients.
Definition at line 183 of file gpuSplineUtils.cuh.
|
private |
GPU arrays to hold TF1 coefficients.
Definition at line 192 of file gpuSplineUtils.cuh.
|
private |
KS: GPU arrays to hold X coefficient.
Definition at line 180 of file gpuSplineUtils.cuh.
|
private |
KS: GPU Number of knots per spline.
Definition at line 186 of file gpuSplineUtils.cuh.
|
private |
KS: GPU map keeping track how many parameters applies to each event, we keep two numbers here {number of splines per event, index where splines start for a given event}.
Definition at line 175 of file gpuSplineUtils.cuh.
|
private |
KS: GPU map keeping track how many parameters applies to each event, we keep two numbers here {number of TF1 per event, index where TF1 start for a given event}.
Definition at line 177 of file gpuSplineUtils.cuh.
|
private |
GPU arrays to hold number of points.
Definition at line 194 of file gpuSplineUtils.cuh.
|
private |
CW: parameter value on GPU.
Definition at line 199 of file gpuSplineUtils.cuh.
|
private |
CW: GPU array with the number of points per spline (not per spline point!)
Definition at line 189 of file gpuSplineUtils.cuh.
|
private |
CW: GPU array with the number of points per TF1 object.
Definition at line 196 of file gpuSplineUtils.cuh.
|
private |
CW: Spline segment on GPU.
Definition at line 201 of file gpuSplineUtils.cuh.
|
private |
GPU arrays to hold weight for event.
Definition at line 204 of file gpuSplineUtils.cuh.
|
private |
GPU arrays to hold weight for each spline.
Definition at line 206 of file gpuSplineUtils.cuh.
|
private |
GPU arrays to hold weight for each TF1.
Definition at line 208 of file gpuSplineUtils.cuh.
|
private |
KS: Textures are L1 cache variables which are well optimised for fetching. Make texture only for variables you often access but rarely overwrite. There are limits on texture memory so don't use huge arrays.
Definition at line 225 of file gpuSplineUtils.cuh.
|
private |
KS: Map keeping track how many parameters applies to each event, we keep two numbers here {number of splines per event, index where splines start for a given event}.
Definition at line 227 of file gpuSplineUtils.cuh.
|
private |
KS: Map keeping track how many parameters applies to each event, we keep two numbers here {number of TF1 per event, index where TF1 start for a given event}.
Definition at line 229 of file gpuSplineUtils.cuh.