MaCh3  2.4.2
Reference Guide
SplineMonolith.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include "Splines/SplineBase.h"
4 
5 //KS: Joy of forward declaration https://gieseanw.wordpress.com/2018/02/25/the-joys-of-forward-declarations-results-from-the-real-world/
6 class SMonolithGPU;
7 
11 class SMonolith : public SplineBase {
12  public:
18  SMonolith(std::vector<std::vector<TResponseFunction_red*> > &MasterSpline,
19  const std::vector<RespFuncType> &SplineType,
20  const bool SaveFlatTree = false,
21  const std::string& _FastSplineName = "SplineFile.root");
24  SMonolith(const std::string& FileName);
26  virtual ~SMonolith();
27 
29  void Evaluate() override;
30 
32  std::string GetName() const override {return "SplineMonolith";};
33 
35  void SynchroniseMemTransfer() const override;
36 
40  const float* retPointer(const int event) const {return &cpu_total_weights[event];}
41 
44  void setSplinePointers(std::vector< const double* > spline_ParsPointers) {
45  splineParsPointer = spline_ParsPointers;
46  for (M3::int_t i = 0; i < nParams; ++i) SplineInfoArray[i].splineParsPointer = spline_ParsPointers[i];
47  };
48 
50  float* cpu_weights;
53 
55  void PrepareSplineFile(std::string FileName) override;
58  void LoadSplineFile(std::string FileName) override;
59 
60  private:
62  inline void Initialise();
73  inline void ScanMasterSpline(std::vector<std::vector<TResponseFunction_red*> > & MasterSpline,
74  unsigned int &nEvents,
75  short int &MaxPoints,
76  short int &numParams,
77  int &nSplines,
78  unsigned int &NSplinesValid,
79  unsigned int &numKnots,
80  unsigned int &nTF1Valid,
81  unsigned int &nTF1_coeff,
82  const std::vector<RespFuncType> &SplineType);
85  inline void PrepareForGPU(std::vector<std::vector<TResponseFunction_red*> > &MasterSpline, const std::vector<RespFuncType> &SplineType);
87  inline void MoveToGPU();
88 
90  inline void PrintInitialsiation();
91 
98  inline void getSplineCoeff_SepMany(TSpline3_red* &spl, int &nPoints, float *&xArray, float *&manyArray);
99 
101  inline void CalcSplineWeights() override;
103  inline void ModifyWeights() override;
105  inline void ModifyWeights_GPU();
106 
108  std::vector< const double* > splineParsPointer;
109 
111  unsigned int NEvents;
113  short int _max_knots;
115  std::vector<int> index_spline_cpu;
117  std::vector<int> index_TF1_cpu;
118 
120  unsigned int NSplines_valid;
122  unsigned int NTF1_valid;
123 
125  unsigned int NSplines_total_large;
126 
128  unsigned int nKnots;
130  unsigned int nTF1coeff;
131 
136 
138  std::vector<unsigned int> cpu_nParamPerEvent;
139 
141  std::vector<unsigned int> cpu_nParamPerEvent_tf1;
142 
145 
148 
150  std::vector<float> cpu_coeff_TF1_many;
151 
153  std::vector<short int> cpu_nPoints_arr;
154 
156  std::vector<short int> cpu_paramNo_TF1_arr;
157 
160 
162  std::string FastSplineName;
163 };
Class responsible for calculating spline weight on GPU.
Even-by-event class calculating response for spline parameters. It is possible to use GPU acceleratio...
SMonolithGPU * gpu_spline_handler
KS: Store info about Spline monolith, this allow to obtain better step time. As all necessary informa...
std::vector< unsigned int > cpu_nParamPerEvent
KS: CPU map keeping track how many parameters applies to each event, we keep two numbers here {number...
unsigned int NTF1_valid
Number of valid TF1.
std::vector< short int > cpu_nPoints_arr
CPU arrays to hold number of points.
void SynchroniseMemTransfer() const override
KS: After calculations are done on GPU we copy memory to CPU. This operation is asynchronous meaning ...
void Evaluate() override
CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; proba...
unsigned int nTF1coeff
Sum of all coefficients over all TF1.
void setSplinePointers(std::vector< const double * > spline_ParsPointers)
KS: Set pointers to spline params.
std::vector< float > cpu_coeff_TF1_many
CPU arrays to hold TF1 coefficients.
void Initialise()
KS: Set everything to null etc.
SplineMonoStruct * cpu_spline_handler
KS: Store info about Spline monolith, this allow to obtain better step time. As all necessary informa...
float * cpu_weights_tf1_var
CPU arrays to hold weight for each TF1.
void PrintInitialsiation()
KS: Print info about how much knots etc has been initialised.
std::vector< unsigned int > cpu_nParamPerEvent_tf1
KS: CPU map keeping track how many parameters applies to each event, we keep two numbers here {number...
std::vector< int > index_TF1_cpu
holds the index for good TF1; don't do unsigned since starts with negative value!
void ModifyWeights_GPU()
Conversion from valid splines to all.
float * cpu_weights
The returned gpu weights, read by the GPU.
float * cpu_weights_spline_var
CPU arrays to hold weight for each spline.
std::vector< int > index_spline_cpu
holds the index for good splines; don't do unsigned since starts with negative value!
void PrepareForGPU(std::vector< std::vector< TResponseFunction_red * > > &MasterSpline, const std::vector< RespFuncType > &SplineType)
CW: Prepare the TSpline3_red objects for the GPU.
void CalcSplineWeights() override
CPU based code which eval weight for each spline.
std::string FastSplineName
Name of Fast Spline to which will be saved.
const float * retPointer(const int event) const
KS: Get pointer to total weight to make fit faster wrooom!
void MoveToGPU()
CW: The shared initialiser from constructors of TResponseFunction_red.
bool SaveSplineFile
Flag telling whether we are saving spline monolith into handy root file.
float * cpu_total_weights
KS: This holds the total CPU weights that gets read in samplePDFND.
unsigned int NSplines_valid
Number of valid splines.
SMonolith(std::vector< std::vector< TResponseFunction_red * > > &MasterSpline, const std::vector< RespFuncType > &SplineType, const bool SaveFlatTree=false, const std::string &_FastSplineName="SplineFile.root")
Constructor.
virtual ~SMonolith()
Destructor for SMonolith class.
void ModifyWeights() override
Calc total event weight.
unsigned int NEvents
Number of events.
void getSplineCoeff_SepMany(TSpline3_red *&spl, int &nPoints, float *&xArray, float *&manyArray)
CW: This loads up coefficients into two arrays: one x array and one yabcd array.
void LoadSplineFile(std::string FileName) override
KS: Load preprocessed spline file.
unsigned int NSplines_total_large
Number of total splines if each event had every parameter's spline.
void ScanMasterSpline(std::vector< std::vector< TResponseFunction_red * > > &MasterSpline, unsigned int &nEvents, short int &MaxPoints, short int &numParams, int &nSplines, unsigned int &NSplinesValid, unsigned int &numKnots, unsigned int &nTF1Valid, unsigned int &nTF1_coeff, const std::vector< RespFuncType > &SplineType)
CW: Function to scan through the MasterSpline of TSpline3.
void PrepareSplineFile(std::string FileName) override
KS: Prepare spline file that can be used for fast loading.
short int _max_knots
Max knots for production.
std::vector< short int > cpu_paramNo_TF1_arr
CW: CPU array with the number of points per spline (not per spline point!)
std::vector< const double * > splineParsPointer
This holds pointer to parameter position which we later copy paste it to GPU.
unsigned int nKnots
Sum of all knots over all splines.
std::string GetName() const override
Get class name.
Base class for calculating weight from spline.
Definition: SplineBase.h:25
short int nParams
Number of parameters that have splines.
Definition: SplineBase.h:81
std::vector< FastSplineInfo > SplineInfoArray
Definition: SplineBase.h:74
CW: Reduced TSpline3 class.
int int_t
Definition: Core.h:38
KS: Struct storing information for spline monolith.
Definition: SplineCommon.h:30