MaCh3  2.5.0
Reference Guide
UnbinnedSplineHandler.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/
7 
12  public:
18  UnbinnedSplineHandler(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  UnbinnedSplineHandler(const std::string& FileName);
26  virtual ~UnbinnedSplineHandler();
27 
29  void Evaluate() final;
30 
32  std::string GetName() const override {return "SplineMonolith";};
33 
35  void SynchroniseMemTransfer() const final;
36 
40  const M3::float_t* RetPointer(const int event) const {return &cpu_total_weights[event];}
41 
44  void setSplinePointers(std::vector< const M3::float_t* > spline_ParsPointers) {
45  for (M3::int_t i = 0; i < nParams; ++i) SplineInfoArray[i].splineParsPointer = spline_ParsPointers[i];
46  };
47 
49  void PrepareSplineFile(std::string FileName) final;
52  void LoadSplineFile(std::string FileName) final;
53  private:
55  void Initialise();
66  void ScanMasterSpline(std::vector<std::vector<TResponseFunction_red*> > & MasterSpline,
67  unsigned int &nEvents,
68  short int &MaxPoints,
69  short int &numParams,
70  int &nSplines,
71  unsigned int &NSplinesValid,
72  unsigned int &numKnots,
73  unsigned int &nTF1Valid,
74  unsigned int &nTF1_coeff,
75  const std::vector<RespFuncType> &SplineType);
78  void PrepareForGPU(std::vector<std::vector<TResponseFunction_red*> > &MasterSpline, const std::vector<RespFuncType> &SplineType);
80  void MoveToGPU();
81  void SetupSegments();
82 
84  void PrintInitialsiation() const;
85 
92  void GetSplineCoeff_SepMany(TSpline3_red* &spl, int &nPoints, float *&xArray, float *&manyArray) const;
93 
95  void CalcSplineWeights() final;
97  void CalcTotalEventWeight();
98 
100  unsigned int NEvents;
102  short int _max_knots;
103 
105  unsigned int NSplines_valid;
107  unsigned int NTF1_valid;
108 
110  unsigned int nKnots;
112  unsigned int nTF1coeff;
113 
120 
122  std::vector<unsigned int> cpu_nParamPerEvent;
123 
125  std::vector<unsigned int> cpu_nParamPerEvent_tf1;
126 
129 
132 
134  std::vector<float> cpu_coeff_TF1_many;
135 
137  std::vector<short int> cpu_nPoints_arr;
138 
140  std::vector<short int> cpu_paramNo_TF1_arr;
141 
144 
146  std::string FastSplineName;
147 };
Base class for calculating weight from spline.
Definition: SplineBase.h:27
short int nParams
Number of parameters that have splines.
Definition: SplineBase.h:81
std::vector< FastSplineInfo > SplineInfoArray
Definition: SplineBase.h:74
Class responsible for calculating spline weight on GPU.
CW: Reduced TSpline3 class.
Even-by-event class calculating response for spline parameters. It is possible to use GPU acceleratio...
float * cpu_weights_tf1_var
CPU arrays to hold weight for each TF1.
float * cpu_weights_spline_var
CPU arrays to hold weight for each spline.
unsigned int nTF1coeff
Sum of all coefficients over all TF1.
unsigned int NTF1_valid
Number of valid TF1.
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) final
KS: Prepare spline file that can be used for fast loading.
std::vector< unsigned int > cpu_nParamPerEvent
KS: CPU map keeping track how many parameters applies to each event, we keep two numbers here {number...
void Evaluate() final
CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; proba...
void SynchroniseMemTransfer() const final
KS: After calculations are done on GPU we copy memory to CPU. This operation is asynchronous meaning ...
const M3::float_t * RetPointer(const int event) const
KS: Get pointer to total weight to make fit faster wrooom!
std::vector< short int > cpu_paramNo_TF1_arr
CW: CPU array with the number of points per spline (not per spline point!)
void PrintInitialsiation() const
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...
SplineMonolithGPU * gpu_spline_handler
KS: Store info about Spline monolith, this allow to obtain better step time. As all necessary informa...
unsigned int NSplines_valid
Number of valid splines.
void MoveToGPU()
CW: The shared initialiser from constructors of TResponseFunction_red.
unsigned int nKnots
Sum of all knots over all splines.
void CalcSplineWeights() final
CPU based code which eval weight for each spline.
void setSplinePointers(std::vector< const M3::float_t * > spline_ParsPointers)
KS: Set pointers to spline params.
short int _max_knots
Max knots for production.
std::vector< short int > cpu_nPoints_arr
CPU arrays to hold number of points.
void PrepareForGPU(std::vector< std::vector< TResponseFunction_red * > > &MasterSpline, const std::vector< RespFuncType > &SplineType)
CW: Prepare the TSpline3_red objects for the GPU.
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...
void LoadSplineFile(std::string FileName) final
KS: Load preprocessed spline file.
virtual ~UnbinnedSplineHandler()
Destructor for UnbinnedSplineHandler class.
std::string GetName() const override
Get class name.
std::vector< float > cpu_coeff_TF1_many
CPU arrays to hold TF1 coefficients.
void GetSplineCoeff_SepMany(TSpline3_red *&spl, int &nPoints, float *&xArray, float *&manyArray) const
CW: This loads up coefficients into two arrays: one x array and one yabcd array.
M3::float_t * cpu_total_weights
KS: This holds the total CPU weights that gets read in SampleHandler.
UnbinnedSplineHandler(std::vector< std::vector< TResponseFunction_red * > > &MasterSpline, const std::vector< RespFuncType > &SplineType, const bool SaveFlatTree=false, const std::string &_FastSplineName="SplineFile.root")
Constructor.
bool SaveSplineFile
Flag telling whether we are saving spline monolith into handy root file.
unsigned int NEvents
Number of events.
std::string FastSplineName
Name of Fast Spline to which will be saved.
void CalcTotalEventWeight()
Calc total event weight.
Main namespace for MaCh3 software.
double float_t
Definition: Core.h:37
int int_t
Definition: Core.h:38
KS: Struct storing information for spline monolith.
Definition: SplineCommon.h:35