MaCh3  2.4.2
Reference Guide
Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | List of all members
SMonolith Class Reference

Even-by-event class calculating response for spline parameters. It is possible to use GPU acceleration. More...

#include <Splines/SplineMonolith.h>

Inheritance diagram for SMonolith:
[legend]
Collaboration diagram for SMonolith:
[legend]

Public Member Functions

 SMonolith (std::vector< std::vector< TResponseFunction_red * > > &MasterSpline, const std::vector< RespFuncType > &SplineType, const bool SaveFlatTree=false, const std::string &_FastSplineName="SplineFile.root")
 Constructor. More...
 
 SMonolith (const std::string &FileName)
 Constructor where you pass path to preprocessed root FileName. More...
 
virtual ~SMonolith ()
 Destructor for SMonolith class. More...
 
void Evaluate () override
 CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; probably the best one here! Same thing but pass parameter spline segments instead of variations. More...
 
std::string GetName () const override
 Get class name. More...
 
void SynchroniseMemTransfer () const override
 KS: After calculations are done on GPU we copy memory to CPU. This operation is asynchronous meaning while memory is being copied some operations are being carried. Memory must be copied before actual reweight. This function make sure all has been copied. More...
 
const float * retPointer (const int event) const
 KS: Get pointer to total weight to make fit faster wrooom! More...
 
void setSplinePointers (std::vector< const double * > spline_ParsPointers)
 KS: Set pointers to spline params. More...
 
void PrepareSplineFile (std::string FileName) override
 KS: Prepare spline file that can be used for fast loading. More...
 
void LoadSplineFile (std::string FileName) override
 KS: Load preprocessed spline file. More...
 
- Public Member Functions inherited from SplineBase
 SplineBase ()
 Constructor. More...
 
virtual ~SplineBase ()
 Destructor. More...
 
short int GetNParams () const
 Get number of spline parameters. More...
 

Public Attributes

float * cpu_weights
 The returned gpu weights, read by the GPU. More...
 
float * cpu_total_weights
 KS: This holds the total CPU weights that gets read in samplePDFND. More...
 

Private Member Functions

void Initialise ()
 KS: Set everything to null etc. More...
 
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. More...
 
void PrepareForGPU (std::vector< std::vector< TResponseFunction_red * > > &MasterSpline, const std::vector< RespFuncType > &SplineType)
 CW: Prepare the TSpline3_red objects for the GPU. More...
 
void MoveToGPU ()
 CW: The shared initialiser from constructors of TResponseFunction_red. More...
 
void PrintInitialsiation ()
 KS: Print info about how much knots etc has been initialised. More...
 
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. More...
 
void CalcSplineWeights () override
 CPU based code which eval weight for each spline. More...
 
void ModifyWeights () override
 Calc total event weight. More...
 
void ModifyWeights_GPU ()
 Conversion from valid splines to all. More...
 

Private Attributes

std::vector< const double * > splineParsPointer
 This holds pointer to parameter position which we later copy paste it to GPU. More...
 
unsigned int NEvents
 Number of events. More...
 
short int _max_knots
 Max knots for production. More...
 
std::vector< int > index_spline_cpu
 holds the index for good splines; don't do unsigned since starts with negative value! More...
 
std::vector< int > index_TF1_cpu
 holds the index for good TF1; don't do unsigned since starts with negative value! More...
 
unsigned int NSplines_valid
 Number of valid splines. More...
 
unsigned int NTF1_valid
 Number of valid TF1. More...
 
unsigned int NSplines_total_large
 Number of total splines if each event had every parameter's spline. More...
 
unsigned int nKnots
 Sum of all knots over all splines. More...
 
unsigned int nTF1coeff
 Sum of all coefficients over all TF1. More...
 
float * cpu_weights_spline_var
 CPU arrays to hold weight for each spline. More...
 
float * cpu_weights_tf1_var
 CPU arrays to hold weight for each TF1. More...
 
std::vector< unsigned int > cpu_nParamPerEvent
 KS: CPU 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...
 
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 of TF1 per event, index where TF1 start for a given event}. More...
 
SplineMonoStructcpu_spline_handler
 KS: Store info about Spline monolith, this allow to obtain better step time. As all necessary information for spline weight calculation are here meaning better cache hits. More...
 
SMonolithGPUgpu_spline_handler
 KS: Store info about Spline monolith, this allow to obtain better step time. As all necessary information for spline weight calculation are here meaning better cache hits. More...
 
std::vector< float > cpu_coeff_TF1_many
 CPU arrays to hold TF1 coefficients. More...
 
std::vector< short int > cpu_nPoints_arr
 CPU arrays to hold number of points. More...
 
std::vector< short int > cpu_paramNo_TF1_arr
 CW: CPU array with the number of points per spline (not per spline point!) More...
 
bool SaveSplineFile
 Flag telling whether we are saving spline monolith into handy root file. More...
 
std::string FastSplineName
 Name of Fast Spline to which will be saved. More...
 

Additional Inherited Members

- Protected Member Functions inherited from SplineBase
void FindSplineSegment ()
 CW:Code used in step by step reweighting, Find Spline Segment for each param. More...
 
void PrepareFastSplineInfoDir (std::unique_ptr< TFile > &SplineFile) const
 KS: Prepare Fast Spline Info within SplineFile. More...
 
void LoadFastSplineInfoDir (std::unique_ptr< TFile > &SplineFile)
 KS: Load preprocessed FastSplineInfo. More...
 
void getTF1Coeff (TF1_red *&spl, int &nPoints, float *&coeffs)
 CW: Gets the polynomial coefficients for TF1. More...
 
- Protected Attributes inherited from SplineBase
std::vector< FastSplineInfoSplineInfoArray
 
short int * SplineSegments
 
float * ParamValues
 Store parameter values they are not in FastSplineInfo as in case of GPU we need to copy paste it to GPU. More...
 
short int nParams
 Number of parameters that have splines. More...
 

Detailed Description

Even-by-event class calculating response for spline parameters. It is possible to use GPU acceleration.

Author
Clarence Wret
Kamil Skwarczynski

Definition at line 11 of file SplineMonolith.h.

Constructor & Destructor Documentation

◆ SMonolith() [1/2]

SMonolith::SMonolith ( std::vector< std::vector< TResponseFunction_red * > > &  MasterSpline,
const std::vector< RespFuncType > &  SplineType,
const bool  SaveFlatTree = false,
const std::string &  _FastSplineName = "SplineFile.root" 
)

Constructor.

Parameters
MasterSplineVector of TSpline3 pointers which we strip back
SplineTypeWhether object is TSpline3 or TF1
SaveFlatTreeWhether we want to save monolith into speedy flat tree
_FastSplineNameName to which spline file will be saved

Definition at line 38 of file SplineMonolith.cpp.

42 : SplineBase() {
43 // *****************************************
44  //KS: If true it will save spline monolith into huge ROOT file
45  SaveSplineFile = SaveFlatTree;
46  FastSplineName = _FastSplineName;
47  Initialise();
48  MACH3LOG_INFO("-- GPUING WITH arrays and master spline containing TResponseFunction_red");
49 
50  // Convert the TSpline3 pointers to the reduced form and call the reduced constructor
51  PrepareForGPU(MasterSpline, SplineType);
52 }
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
void Initialise()
KS: Set everything to null etc.
void PrepareForGPU(std::vector< std::vector< TResponseFunction_red * > > &MasterSpline, const std::vector< RespFuncType > &SplineType)
CW: Prepare the TSpline3_red objects for the GPU.
std::string FastSplineName
Name of Fast Spline to which will be saved.
bool SaveSplineFile
Flag telling whether we are saving spline monolith into handy root file.
SplineBase()
Constructor.
Definition: SplineBase.cpp:7

◆ SMonolith() [2/2]

SMonolith::SMonolith ( const std::string &  FileName)

Constructor where you pass path to preprocessed root FileName.

Parameters
FileNamepath to pre-processed root file containing stripped monolith info

Definition at line 502 of file SplineMonolith.cpp.

503  : SplineBase() {
504 // *****************************************
505  Initialise();
506  MACH3LOG_INFO("-- GPUING WITH {X} and {Y,B,C,D} arrays and master spline containing TSpline3_red");
507  // Convert the TSpline3 pointers to the reduced form and call the reduced constructor
508  LoadSplineFile(FileName);
509 }
void LoadSplineFile(std::string FileName) override
KS: Load preprocessed spline file.

◆ ~SMonolith()

SMonolith::~SMonolith ( )
virtual

Destructor for SMonolith class.

Definition at line 683 of file SplineMonolith.cpp.

683  {
684 // *****************************************
685  #ifdef MaCh3_CUDA
687  #ifndef Weight_On_SplineBySpline_Basis
689  #endif
690  );
691 
692  //KS: Since we declared them using CUDA alloc we have to free memory using also cuda functions
694 
695  delete gpu_spline_handler;
696  #else
697  if(SplineSegments != nullptr) delete[] SplineSegments;
698  if(ParamValues != nullptr) delete[] ParamValues;
699  if(cpu_total_weights != nullptr) delete[] cpu_total_weights;
700  #endif
701 
702  if(cpu_weights != nullptr) delete[] cpu_weights;
703  if(cpu_weights_spline_var != nullptr) delete[] cpu_weights_spline_var;
704  if(cpu_weights_tf1_var != nullptr) delete[] cpu_weights_tf1_var;
705 
706  if(cpu_spline_handler != nullptr) delete cpu_spline_handler;
707 }
__host__ void CleanupGPU_SplineMonolith(float *cpu_total_weights)
This function deallocates the resources allocated for the separate {x} and {ybcd} arrays in the and T...
__host__ void CleanupGPU_Segments(short int *segment, float *vals)
Clean up pinned variables at CPU.
SMonolithGPU * gpu_spline_handler
KS: Store info about Spline monolith, this allow to obtain better step time. As all necessary informa...
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.
float * cpu_weights
The returned gpu weights, read by the GPU.
float * cpu_weights_spline_var
CPU arrays to hold weight for each spline.
float * cpu_total_weights
KS: This holds the total CPU weights that gets read in samplePDFND.
short int * SplineSegments
Definition: SplineBase.h:77
float * ParamValues
Store parameter values they are not in FastSplineInfo as in case of GPU we need to copy paste it to G...
Definition: SplineBase.h:79

Member Function Documentation

◆ CalcSplineWeights()

void SMonolith::CalcSplineWeights ( )
inlineoverrideprivatevirtual

CPU based code which eval weight for each spline.

Implements SplineBase.

Definition at line 805 of file SplineMonolith.cpp.

805  {
806 //*********************************************************
807  #ifdef MULTITHREAD
808  //KS: Open parallel region
809  #pragma omp parallel
810  {
811  #endif
812  //KS: First we calculate
813  #ifdef MULTITHREAD
814  #pragma omp for simd nowait
815  #endif
816  for (unsigned int splineNum = 0; splineNum < NSplines_valid; ++splineNum)
817  {
818  //CW: Which Parameter we are accessing
819  const short int Param = cpu_spline_handler->paramNo_arr[splineNum];
820 
821  //CW: Avoids doing costly binary search on GPU
822  const short int segment = SplineSegments[Param];
823 
824  //KS: Segment for coeff_x is simply parameter*max knots + segment as each parameters has the same spacing
825  const short int segment_X = short(Param*_max_knots+segment);
826 
827  //KS: Find knot position in out monolithical structure
828  const unsigned int CurrentKnotPos = cpu_spline_handler->nKnots_arr[splineNum]*_nCoeff_+segment*_nCoeff_;
829 
830  // We've read the segment straight from CPU and is saved in segment_gpu
831  // polynomial parameters from the monolithic splineMonolith
832  const float fY = cpu_spline_handler->coeff_many[CurrentKnotPos];
833  const float fB = cpu_spline_handler->coeff_many[CurrentKnotPos + 1];
834  const float fC = cpu_spline_handler->coeff_many[CurrentKnotPos + 2];
835  const float fD = cpu_spline_handler->coeff_many[CurrentKnotPos + 3];
836  // The is the variation itself (needed to evaluate variation - stored spline point = dx)
837  const float dx = ParamValues[Param] - cpu_spline_handler->coeff_x[segment_X];
838 
839  //CW: Wooow, let's use some fancy intrinsic and pull down the processing time by <1% from normal multiplication! HURRAY
840  cpu_weights_spline_var[splineNum] = fmaf(dx, fmaf(dx, fmaf(dx, fD, fC), fB), fY);
841  // Or for the more "easy to read" version:
842  //cpu_weights_spline_var[splineNum] = (fY+dx*(fB+dx*(fC+dx*fD)));
843  }
844 
845  #ifdef MULTITHREAD
846  #pragma omp for simd
847  #endif
848  for (unsigned int tf1Num = 0; tf1Num < NTF1_valid; ++tf1Num)
849  {
850  // The is the variation itself (needed to evaluate variation - stored spline point = dx)
851  const float x = ParamValues[cpu_paramNo_TF1_arr[tf1Num]];
852 
853  // Read the coefficients
854  const unsigned int TF1_Index = tf1Num * _nTF1Coeff_;
855  const float a = cpu_coeff_TF1_many[TF1_Index];
856  const float b = cpu_coeff_TF1_many[TF1_Index + 1];
857 
858  cpu_weights_tf1_var[tf1Num] = fmaf(a, x, b);
859  // cpu_weights_tf1_var[tf1Num] = a*x + b;
860  //cpu_weights_tf1_var[splineNum] = 1 + a*x + b*x*x + c*x*x*x + d*x*x*x*x + e*x*x*x*x*x;
861  }
862  #ifdef MULTITHREAD
863  //KS: End parallel region
864  }
865  #endif
866 }
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...
Definition: SplineCommon.h:13
constexpr int _nTF1Coeff_
KS: For TF1 we store at most 5 coefficients, we could make it more flexible but for now define it her...
Definition: SplineCommon.h:15
unsigned int NTF1_valid
Number of valid TF1.
std::vector< float > cpu_coeff_TF1_many
CPU arrays to hold TF1 coefficients.
unsigned int NSplines_valid
Number of valid splines.
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< unsigned int > nKnots_arr
KS: CPU Number of knots per spline.
Definition: SplineCommon.h:41
std::vector< float > coeff_x
KS: CPU arrays to hold X coefficient.
Definition: SplineCommon.h:32
std::vector< float > coeff_many
CPU arrays to hold other coefficients.
Definition: SplineCommon.h:38
std::vector< short int > paramNo_arr
CW: CPU array with the number of points per spline (not per spline point!)
Definition: SplineCommon.h:44

◆ Evaluate()

void SMonolith::Evaluate ( )
overridevirtual

CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; probably the best one here! Same thing but pass parameter spline segments instead of variations.

Implements SplineBase.

Definition at line 790 of file SplineMonolith.cpp.

790  {
791 // *****************************************
792  // There's a parameter mapping that goes from spline parameter to a global parameter index
793  // Find the spline segments
795 
796  //KS: Huge MP loop over all valid splines
798 
799  //KS: Huge MP loop over all events calculating total weight
800  ModifyWeights();
801 }
void CalcSplineWeights() override
CPU based code which eval weight for each spline.
void ModifyWeights() override
Calc total event weight.
void FindSplineSegment()
CW:Code used in step by step reweighting, Find Spline Segment for each param.
Definition: SplineBase.cpp:24

◆ GetName()

std::string SMonolith::GetName ( ) const
inlineoverridevirtual

Get class name.

Reimplemented from SplineBase.

Definition at line 32 of file SplineMonolith.h.

32 {return "SplineMonolith";};

◆ getSplineCoeff_SepMany()

void SMonolith::getSplineCoeff_SepMany ( TSpline3_red *&  spl,
int &  nPoints,
float *&  xArray,
float *&  manyArray 
)
inlineprivate

CW: This loads up coefficients into two arrays: one x array and one yabcd array.

CW: This should maximize our cache hits!

Parameters
splpointer to TSpline3_red
nPointsnumber of knots
xArrayarray X value for each knot
manyArrayArray holding coefficients for each knot

Definition at line 713 of file SplineMonolith.cpp.

713  {
714 // *****************************************
715  // Initialise all arrays to 1.0
716  for (int i = 0; i < _max_knots; ++i) {
717  xArray[i] = 1.0;
718  for (int j = 0; j < _nCoeff_; j++) {
719  manyArray[i*_nCoeff_+j] = 1.0;
720  }
721  }
722  // Get number of points in spline
723  int Np = spl->GetNp();
724  // If spline is flat, set number of knots to 1.0,
725  // This is used later to expedite the calculations for flat splines
726  // tmpArray[0] is number of knots
727  nPoints = Np;
728  if (Np > _max_knots) {
729  MACH3LOG_ERROR("Error, number of points is greater than saved {}", _max_knots);
730  MACH3LOG_ERROR("This _WILL_ cause problems with GPU splines and _SHOULD_ be fixed!");
731  MACH3LOG_ERROR("nPoints = {}, _max_knots = {}", nPoints, _max_knots);
732  throw MaCh3Exception(__FILE__ , __LINE__ );
733  }
734 
735  // The coefficients we're writing to
736  M3::float_t x, y, b, c, d;
737  // TSpline3 can only take doubles, not floats
738  // But our GPU is slow with doubles, so need to cast to float
739  for(int i = 0; i < Np; i++) {
740  // Get the coefficients from the TSpline3 object
741  spl->GetCoeff(i, x, y, b, c, d);
742  // Write the arrays
743  xArray[i] = float(x);
744  manyArray[i*_nCoeff_] = float(y); // 4 because manyArray stores y,b,c,d
745  manyArray[i*_nCoeff_+1] = float(b);
746  manyArray[i*_nCoeff_+2] = float(c);
747  manyArray[i*_nCoeff_+3] = float(d);
748  if((xArray[i] == -999) || (manyArray[i*_nCoeff_] == -999) || (manyArray[i*_nCoeff_ +1] == -999) || (manyArray[i*_nCoeff_+2] == -999) || (manyArray[i*_nCoeff_+3] == -999)){
749  MACH3LOG_ERROR("*********** Bad params in getSplineCoeff_SepMany() ************");
750  MACH3LOG_ERROR("pre cast to float (x, y, b, c, d) = {:.2f}, {:.2f}, {:.2f}, {:.2f}, {:.2f}", x, y, b, c, d);
751  MACH3LOG_ERROR("pre cast to float (x, y, b, c, d) = {:.2f}, {:.2f}, {:.2f}, {:.2f}, {:.2f}", xArray[i], manyArray[i*_nCoeff_], manyArray[i*_nCoeff_+1], manyArray[i*_nCoeff_+2], manyArray[i*_nCoeff_+3]);
752  MACH3LOG_ERROR("This will cause problems when preparing for GPU");
753  MACH3LOG_ERROR("***************************************************************");
754  }
755  }
756 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
Custom exception class used throughout MaCh3.
M3::int_t GetNp() override
CW: Get the number of points.
void GetCoeff(int segment, M3::float_t &x, M3::float_t &y, M3::float_t &b, M3::float_t &c, M3::float_t &d)
CW: Get the coefficient of a given segment.
double float_t
Definition: Core.h:37

◆ Initialise()

void SMonolith::Initialise ( )
inlineprivate

KS: Set everything to null etc.

Definition at line 12 of file SplineMonolith.cpp.

12  {
13 // *****************************************
14 #ifdef MaCh3_CUDA
15  MACH3LOG_INFO("Using GPU version event by event monolith");
16  gpu_spline_handler = nullptr;
17 #endif
18 
20 
21  nKnots = 0;
22  nTF1coeff = 0;
23  NEvents = 0;
24  _max_knots = 0;
25 
26  NSplines_valid = 0;
27  NTF1_valid = 0;
29 
30  cpu_weights_spline_var = nullptr;
31  cpu_weights = nullptr;
32  cpu_weights_tf1_var = nullptr;
33 
34  cpu_total_weights = nullptr;
35 }
unsigned int nTF1coeff
Sum of all coefficients over all TF1.
unsigned int NEvents
Number of events.
unsigned int NSplines_total_large
Number of total splines if each event had every parameter's spline.
unsigned int nKnots
Sum of all knots over all splines.
KS: Struct storing information for spline monolith.
Definition: SplineCommon.h:30

◆ LoadSplineFile()

void SMonolith::LoadSplineFile ( std::string  FileName)
overridevirtual

KS: Load preprocessed spline file.

Parameters
FileNamePath to ROOT file with predefined reduced Spline Monolith

Implements SplineBase.

Definition at line 513 of file SplineMonolith.cpp.

513  {
514 // *****************************************
515  #ifdef Weight_On_SplineBySpline_Basis
516  MACH3LOG_ERROR("Trying to load Monolith from file using weight by weight base, this is not supported right now, sorry");
517  throw MaCh3Exception(__FILE__ , __LINE__ );
518  #endif
519 
520  M3::AddPath(FileName);
521  auto SplineFile = std::make_unique<TFile>(FileName.c_str(), "OPEN");
522  TTree *Settings = SplineFile->Get<TTree>("Settings");
523  TTree *Monolith_TF1 = SplineFile->Get<TTree>("Monolith_TF1");
524  TTree *EventInfo = SplineFile->Get<TTree>("EventInfo");
525  TTree *SplineTree = SplineFile->Get<TTree>("SplineTree");
526 
527  unsigned int NEvents_temp;
528  short int nParams_temp;
529  int _max_knots_temp;
530  unsigned int nKnots_temp;
531  unsigned int NSplines_valid_temp;
532  unsigned int nTF1Valid_temp;
533  unsigned int nTF1coeff_temp;
534 
535  Settings->SetBranchAddress("NEvents", &NEvents_temp);
536  Settings->SetBranchAddress("nParams", &nParams_temp);
537  Settings->SetBranchAddress("_max_knots", &_max_knots_temp);
538  Settings->SetBranchAddress("nKnots", &nKnots_temp);
539  Settings->SetBranchAddress("NSplines_valid", &NSplines_valid_temp);
540  Settings->SetBranchAddress("NTF1_valid", &nTF1Valid_temp);
541  Settings->SetBranchAddress("nTF1coeff", &nTF1coeff_temp);
542 
543  Settings->GetEntry(0);
544 
545  NEvents = NEvents_temp;
546  nParams = nParams_temp;
547  _max_knots = static_cast<short int>(_max_knots_temp);
548  nKnots = nKnots_temp;
549  NSplines_valid = NSplines_valid_temp;
550  NTF1_valid = nTF1Valid_temp;
551  nTF1coeff = nTF1coeff_temp;
552 
553  //KS: Since we are going to copy it each step use fancy CUDA memory allocation
554 #ifdef MaCh3_CUDA
557 #else
558  SplineSegments = new short int[nParams]();
559  ParamValues = new float[nParams]();
560 #endif
561 
562  cpu_nParamPerEvent.resize(2*NEvents);
565 
566  //KS: This is tricky as this variable use both by CPU and GPU, however if use CUDA we use cudaMallocHost
567 #ifndef MaCh3_CUDA
568  cpu_total_weights = new float[NEvents]();
569  cpu_weights_spline_var = new float[NSplines_valid]();
570  cpu_weights_tf1_var = new float[NTF1_valid]();
571 #endif
572 
573  SplineTree->SetBranchAddress("SplineObject", &cpu_spline_handler);
574  SplineTree->GetEntry(0);
575 
576  float coeff_tf1 = 0.;
577  Monolith_TF1->SetBranchAddress("cpu_coeff_TF1_many", &coeff_tf1);
578  for(unsigned int i = 0; i < nTF1coeff; i++)
579  {
580  Monolith_TF1->GetEntry(i);
581  cpu_coeff_TF1_many[i] = coeff_tf1;
582  }
583 
584  unsigned int nParamPerEvent = 0;
585  unsigned int nParamPerEvent_tf1 = 0;
586 
587  EventInfo->SetBranchAddress("cpu_nParamPerEvent", &nParamPerEvent);
588  EventInfo->SetBranchAddress("cpu_nParamPerEvent_tf1", &nParamPerEvent_tf1);
589  for(unsigned int i = 0; i < 2*NEvents; i++)
590  {
591  EventInfo->GetEntry(i);
592  cpu_nParamPerEvent[i] = nParamPerEvent;
593  cpu_nParamPerEvent_tf1[i] = nParamPerEvent_tf1;
594  }
595 
596  LoadFastSplineInfoDir(SplineFile);
597 
598  SplineFile->Close();
599 
600  // Print some info; could probably make this to a separate function
602 
603  MoveToGPU();
604 }
__host__ void InitGPU_Vals(float **vals)
Allocate memory for spline segments.
__host__ void InitGPU_Segments(short int **segment)
Allocate memory for spline segments.
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 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...
void MoveToGPU()
CW: The shared initialiser from constructors of TResponseFunction_red.
short int nParams
Number of parameters that have splines.
Definition: SplineBase.h:81
void LoadFastSplineInfoDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed FastSplineInfo.
Definition: SplineBase.cpp:159
void AddPath(std::string &FilePath)
Prepends the MACH3 environment path to FilePath if it is not already present.
Definition: Monitor.cpp:382
Stores info about each MC event used during reweighting routine.

◆ ModifyWeights()

void SMonolith::ModifyWeights ( )
inlineoverrideprivatevirtual

Calc total event weight.

Implements SplineBase.

Definition at line 870 of file SplineMonolith.cpp.

870  {
871 //*********************************************************
872 #ifndef Weight_On_SplineBySpline_Basis
873  #ifdef MULTITHREAD
874  #pragma omp parallel for
875  #endif
876  for (unsigned int EventNum = 0; EventNum < NEvents; ++EventNum)
877  {
878  float totalWeight = 1.0f; // Initialize total weight for each event
879 
880  const unsigned int Offset = 2 * EventNum;
881 
882  // Extract the parameters for the current event
883  const unsigned int startIndex = cpu_nParamPerEvent[Offset + 1];
884  const unsigned int numParams = cpu_nParamPerEvent[Offset];
885 
886  // Compute total weight for the current event
887  #ifdef MULTITHREAD
888  #pragma omp simd reduction(*:totalWeight)
889  #endif
890  for (unsigned int id = 0; id < numParams; ++id) {
891  totalWeight *= cpu_weights_spline_var[startIndex + id];
892  }
893  //Now TF1
894  // Extract the parameters for the current event
895  const unsigned int startIndex_tf1 = cpu_nParamPerEvent_tf1[Offset + 1];
896  const unsigned int numParams_tf1 = cpu_nParamPerEvent_tf1[Offset];
897 
898  // Compute total weight for the current event
899  #ifdef MULTITHREAD
900  #pragma omp simd reduction(*:totalWeight)
901  #endif
902  for (unsigned int id = 0; id < numParams_tf1; ++id) {
903  totalWeight *= cpu_weights_tf1_var[startIndex_tf1 + id];
904  }
905 
906  // Store the total weight for the current event
907  cpu_total_weights[EventNum] = totalWeight;
908  }
909 #else
910  //KS: Name is confusing but what it does it make a nice mapping used for debugging
912 #endif
913 }
void ModifyWeights_GPU()
Conversion from valid splines to all.

◆ ModifyWeights_GPU()

void SMonolith::ModifyWeights_GPU ( )
inlineprivate

Conversion from valid splines to all.

Definition at line 917 of file SplineMonolith.cpp.

917  {
918 //*********************************************************
919 #ifdef Weight_On_SplineBySpline_Basis
920  // Multi-thread here because _numIndex is really quite large!
921  #ifdef MULTITHREAD
922  #pragma omp parallel for
923  #endif
924  for (unsigned int i = 0; i < NSplines_total_large; ++i) {
925  if (index_spline_cpu[i] >= 0) {
927  } else if (index_TF1_cpu[i] >= 0) {
929  } else {
930  cpu_weights[i] = 1.;
931  }
932  }
933 #endif
934 }
std::vector< int > index_TF1_cpu
holds the index for good TF1; don't do unsigned since starts with negative value!
std::vector< int > index_spline_cpu
holds the index for good splines; don't do unsigned since starts with negative value!

◆ MoveToGPU()

void SMonolith::MoveToGPU ( )
inlineprivate

CW: The shared initialiser from constructors of TResponseFunction_red.

Definition at line 305 of file SplineMonolith.cpp.

305  {
306 // *****************************************
307  #ifdef MaCh3_CUDA
308  unsigned int event_size_max = _max_knots * nParams;
309  MACH3LOG_INFO("Total size = {:.2f} MB memory on CPU to move to GPU",
310  (double(sizeof(float) * nKnots * _nCoeff_) + double(sizeof(float) * event_size_max) / 1.E6 +
311  double(sizeof(short int) * NSplines_valid)) / 1.E6);
312  MACH3LOG_INFO("Total TF1 size = {:.2f} MB memory on CPU to move to GPU",
313  double(sizeof(float) * NTF1_valid * _nTF1Coeff_) / 1.E6);
314  MACH3LOG_INFO("GPU weight array (GPU->CPU every step) = {:.2f} MB", static_cast<double>(sizeof(float)) * (NSplines_valid + NTF1_valid) / 1.0e6);
315  #ifndef Weight_On_SplineBySpline_Basis
316  MACH3LOG_INFO("Since you are running Total event weight mode then GPU weight array (GPU->CPU every step) = {:.2f} MB",
317  double(sizeof(float) * NEvents) / 1.E6);
318  #endif
319  MACH3LOG_INFO("Parameter value array (CPU->GPU every step) = {:.4f} MB", double(sizeof(float) * nParams) / 1.E6);
320  //CW: With the new set-up we have: 1 coefficient array of size coeff_array_size, all same size
321  // 1 coefficient array of size coeff_array_size*4, holding y,b,c,d in order (y11,b11,c11,d11; y12,b12,c12,d12;...) where ynm is n = spline number, m = spline point. Should really make array so that order is (y11,b11,c11,d11; y21,b21,c21,d21;...) because it will optimise cache hits I think; try this if you have time
322  // return gpu_weights
323 
325 
326  // The gpu_XY arrays don't actually need initialising, since they are only placeholders for what we'll move onto the GPU. As long as we cudaMalloc the size of the arrays correctly there shouldn't be any problems
327  // Can probably make this a bit prettier but will do for now
328  // Could be a lot smaller of a function...
330  #ifndef Weight_On_SplineBySpline_Basis
332  NEvents,
333  #endif
334  nKnots, // How many entries in coefficient array (*4 for the "many" array)
335  NSplines_valid, // What's the number of splines we have (also number of entries in gpu_nPoints_arr)
336  NTF1_valid,
337  event_size_max //Knots times event number of unique splines
338  );
339 
340  // Move number of splines and spline size to constant GPU memory; every thread does not need a copy...
341  // The implementation lives in splines/gpuSplineUtils.cu
342  // The GPU splines don't actually need declaring but is good for demonstration, kind of
343  // fixed by passing const reference
346 
347  // TFI related now
350  #ifndef Weight_On_SplineBySpline_Basis
351  NEvents,
354  #endif
355  nParams,
357  _max_knots,
358  nKnots,
359  NTF1_valid);
360 
361  // Delete all the coefficient arrays from the CPU once they are on the GPU
364  #ifndef Weight_On_SplineBySpline_Basis
367  #endif
368  delete cpu_spline_handler;
369  cpu_spline_handler = nullptr;
370  MACH3LOG_INFO("Good GPU loading");
371  #endif
372 }
void CleanVector(T &)
Base case: do nothing for non-vector types.
Class responsible for calculating spline weight on GPU.
__host__ void InitGPU_SplineMonolith(float **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.
__host__ void CopyToGPU_SplineMonolith(SplineMonoStruct *cpu_spline_handler, std::vector< float > cpu_many_array_TF1, std::vector< short int > cpu_paramNo_arr_TF1, int n_events, std::vector< unsigned int > cpu_nParamPerEvent, std::vector< unsigned int > cpu_nParamPerEvent_TF1, int n_params, unsigned int n_splines, short int spline_size, unsigned int total_nknots, unsigned int n_tf1)
Copies data from CPU to GPU for the spline monolith.

◆ PrepareForGPU()

void SMonolith::PrepareForGPU ( std::vector< std::vector< TResponseFunction_red * > > &  MasterSpline,
const std::vector< RespFuncType > &  SplineType 
)
inlineprivate

CW: Prepare the TSpline3_red objects for the GPU.

Parameters
MasterSplineVector of TResponseFunction_red pointers which we strip back

Definition at line 56 of file SplineMonolith.cpp.

56  {
57 // *****************************************
58 
59  // Scan for the max number of knots, the number of events (number of splines), and number of parameters
60  int maxnSplines = 0;
61  ScanMasterSpline(MasterSpline,
62  NEvents,
63  _max_knots,
64  nParams,
65  maxnSplines,
67  nKnots,
68  NTF1_valid,
69  nTF1coeff,
70  SplineType);
71 
72  MACH3LOG_INFO("Found {} events", NEvents);
73  MACH3LOG_INFO("Found {} knots at max", _max_knots);
74  MACH3LOG_INFO("Found {} parameters", nParams);
75  MACH3LOG_INFO("Found {} maximum number of splines in an event", maxnSplines);
76  MACH3LOG_INFO("Found total {} knots in all splines", nKnots);
77  MACH3LOG_INFO("Number of splines = {}", NSplines_valid);
78  MACH3LOG_INFO("Found total {} coeffs in all TF1", nTF1coeff);
79  MACH3LOG_INFO("Number of TF1 = {}", NTF1_valid);
80 
81  // Can pass the spline segments to the GPU instead of the values
82  // Make these here and only refill them for each loop, avoiding unnecessary new/delete on each reconfigure
83  //KS: Since we are going to copy it each step use fancy CUDA memory allocation
84  #ifdef MaCh3_CUDA
87  #else
88  SplineSegments = new short int[nParams];
89  ParamValues = new float[nParams];
90  #endif
91 
92  for (M3::int_t j = 0; j < nParams; j++)
93  {
94  SplineSegments[j] = 0;
95  ParamValues[j] = -999;
96  }
97 
98  // Number of objects we have in total if each event has *EVERY* spline. Needed for some arrays
100 
101  unsigned int event_size_max = _max_knots * nParams;
102  // Declare the {x}, {y,b,c,d} arrays for all possible splines which the event has
103  // We'll filter off the flat and "disabled" (e.g. CCQE event should not have MARES spline) ones in the next for loop, but need to declare these beasts here
104 
105  // Declare the {y,b,c,d} for each knot
106  // float because GPU precision (could change to double, but will incur significant speed reduction on GPU unless you're very rich!)
107  cpu_spline_handler->coeff_many.resize(nKnots*_nCoeff_); // *4 because we store y,b,c,d parameters in this array
108  //KS: For x coeff we assume that for given dial (MAQE) spacing is identical, here we are sloppy and assume each dial has the same number of knots, not a big problem
109  cpu_spline_handler->coeff_x.resize(event_size_max);
110 
111  // Set all the big arrays to -999 to keep us safe...
112  for (unsigned int j = 0; j < event_size_max; j++) {
113  cpu_spline_handler->coeff_x[j] = -999;
114  }
115 
116  //CW: With TF1 we only save the coefficients and the order of the polynomial
117  // Makes most sense to have one large monolithic array, but then it becomes impossible to tell apart a coefficient from a "number of points". So have two arrays: one of coefficients and one of number of points
118  // Let's first assume all are of _max_knots size
119  // Now declare the arrays for each point in the valid splines which the event actually has (i.e. include the splines that the event undergoes)
120  // Also make array with the number of points per spline (not per spline point!)
121  // float because GPU precision (could change to double, but will incur significant speed reduction on GPU unless you're very rich!)
122  cpu_nPoints_arr.resize(NTF1_valid);
123  cpu_coeff_TF1_many.resize(nTF1coeff); // *5 because this array holds a,b,c,d,e parameters
124 
125  #ifdef Weight_On_SplineBySpline_Basis
126  // This holds the index of each spline
129 
130  #ifdef MULTITHREAD
131  #pragma omp parallel for
132  #endif
133  for (unsigned int j = 0; j < NSplines_total_large; j++) {
134  index_spline_cpu[j] = -1;
135  index_TF1_cpu[j] = -1;
136  }
137  // This holds the total CPU weights that gets read in samplePDFND
138  cpu_weights = new float[NSplines_total_large];
139  #else
140  //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}
141  cpu_nParamPerEvent.resize(2*NEvents);
143  #ifdef MULTITHREAD
144  #pragma omp parallel for
145  #endif
146  for (unsigned int j = 0; j < 2*NEvents; j++) {
147  cpu_nParamPerEvent[j] = -1;
148  cpu_nParamPerEvent_tf1[j] = -1;
149  }
150  #endif
151 
152  // Make array with the number of points per spline (not per spline point!)
154  //KS: And array which tells where each spline stars in a big monolith array, sort of knot map
157 
158  // Temporary arrays to hold the coefficients for each spline
159  // We get one x, one y, one b,... for each point, so only need to be _max_knots big
160  //KS: Some params has less splines but this is all right main array will get proper number while this temp will be deleted
161  float *x_tmp = new float[_max_knots]();
162  float *many_tmp = new float[_max_knots*_nCoeff_]();
163  float *temp_coeffs = new float[_nTF1Coeff_]();
164 
165  // Count the number of events
166  unsigned int KnotCounter = 0;
167  unsigned int TF1PointsCounter = 0;
168  unsigned int NSplinesCounter = 0;
169  unsigned int TF1sCounter = 0;
170  int ParamCounter = 0;
171  int ParamCounterGlobal = 0;
172  int ParamCounter_TF1 = 0;
173  int ParamCounterGlobalTF1 = 0;
174  // Loop over events and extract the spline coefficients
175  for(unsigned int EventCounter = 0; EventCounter < MasterSpline.size(); ++EventCounter) {
176  // Structure of MasterSpline is std::vector<std::vector<TSpline3*>>
177  // A conventional iterator to count which parameter a given spline should be applied to
178  for(unsigned int ParamNumber = 0; ParamNumber < MasterSpline[EventCounter].size(); ++ParamNumber) {
179  // If NULL we don't have this spline for the event, so move to next spline
180  if (MasterSpline[EventCounter][ParamNumber] == NULL) continue;
181 
182  if(SplineType[ParamNumber] == kTSpline3_red)
183  {
184  //KS: how much knots each spline has
185  int nPoints_tmp = 0;
186  // Get a pointer to the current spline for this event
187  TResponseFunction_red* TespFunc = MasterSpline[EventCounter][ParamNumber];
188  TSpline3_red* CurrSpline = static_cast<TSpline3_red*>(TespFunc);
189 
190  // If the number of knots are greater than 2 the spline is not a dummy and we should extract coefficients to load onto the GPU
191  getSplineCoeff_SepMany(CurrSpline, nPoints_tmp, x_tmp, many_tmp);
192 
193  //KS: One knot means flat spline so ignore
194  if (nPoints_tmp == 1) continue;
195  for (int j = 0; j < _max_knots; ++j) {
196  cpu_spline_handler->coeff_x[ParamNumber*_max_knots + j] = x_tmp[j];
197  }
198  //KS: Contrary to X coeff we keep for other coeff only filled knots, there is no much gain for doing so for x coeff
199  for (int j = 0; j < nPoints_tmp; ++j) {
200  for (int k = 0; k < _nCoeff_; k++) {
201  cpu_spline_handler->coeff_many[KnotCounter*_nCoeff_ + j*_nCoeff_ + k] = many_tmp[j*_nCoeff_+k];
202  }
203  }
204  // Set the parameter number for this spline
205  cpu_spline_handler->paramNo_arr[NSplinesCounter] = short(ParamNumber);
206  //KS: Fill map when each spline starts
207  cpu_spline_handler->nKnots_arr[NSplinesCounter] = KnotCounter;
208  KnotCounter += nPoints_tmp;
209 
210  #ifdef Weight_On_SplineBySpline_Basis
211  // Set the index of the spline so we can tell apart from flat splines
212  index_spline_cpu[EventCounter*nParams + ParamNumber] = NSplinesCounter;
213  #else
214  ++ParamCounter;
215  #endif
216  // Increment the counter for the number of good splines we have
217  ++NSplinesCounter;
218  }
219  else if (SplineType[ParamNumber] == kTF1_red)
220  {
221  // Don't actually use this ever -- we give each spline the maximum number of points found in all splines
222  int nPoints_tmp = 0;
223  // Get a pointer to the current spline for this event
224  TF1_red* CurrSpline = dynamic_cast<TF1_red*>(MasterSpline[EventCounter][ParamNumber]);
225 
226  // If the number of knots are greater than 2 the spline is not a dummy and we should extract coefficients to load onto the GPU
227  getTF1Coeff(CurrSpline, nPoints_tmp, temp_coeffs);
228  for (int j = 0; j < _nTF1Coeff_; ++j) {
229  cpu_coeff_TF1_many[TF1PointsCounter+j] = temp_coeffs[j];
230  }
231  // Save the number of points for this spline
232  cpu_nPoints_arr[TF1sCounter] = short(nPoints_tmp);
233 
234  TF1PointsCounter += nPoints_tmp;
235  // Set the parameter number for this spline
236  cpu_paramNo_TF1_arr[TF1sCounter] = short(ParamNumber);
237  #ifdef Weight_On_SplineBySpline_Basis
238  // Set the index of the spline so we can tell apart from flat splines
239  index_TF1_cpu[EventCounter*nParams + ParamNumber] = TF1sCounter;
240  #else
241  ++ParamCounter_TF1;
242  #endif
243  // Increment the counter for the number of good splines we have
244  ++TF1sCounter;
245  }
246  //KS: Don't delete in debug
247  #ifndef DEBUG
248  delete MasterSpline[EventCounter][ParamNumber];
249  MasterSpline[EventCounter][ParamNumber] = NULL;
250  #endif
251  } // End the loop over the parameters in the MasterSpline
252  #ifndef Weight_On_SplineBySpline_Basis
253  cpu_nParamPerEvent[2*EventCounter] = ParamCounter;
254  cpu_nParamPerEvent[2*EventCounter+1] = ParamCounterGlobal;
255  ParamCounterGlobal += ParamCounter;
256 
257  cpu_nParamPerEvent_tf1[2*EventCounter] = ParamCounter_TF1;
258  cpu_nParamPerEvent_tf1[2*EventCounter+1] = ParamCounterGlobalTF1;
259  ParamCounterGlobalTF1 += ParamCounter_TF1;
260 
261  ParamCounter = 0;
262  ParamCounter_TF1 = 0;
263  #endif
264  } // End the loop over the number of events
265  delete[] many_tmp;
266  delete[] x_tmp;
267  delete[] temp_coeffs;
268 
269  int BadXCounter = 0;
270  for (unsigned int j = 0; j < event_size_max; j++) {
271  if (cpu_spline_handler->coeff_x[j] == -999) BadXCounter++;
272  // Perform checks that all entries have been modified from initial values
273  if (cpu_spline_handler->coeff_x[j] == -999 && BadXCounter < 5) {
274  MACH3LOG_WARN("***** BAD X !! *****");
275  MACH3LOG_WARN("Indicates some parameter doesn't have a single spline");
276  MACH3LOG_WARN("j = {}", j);
277  //throw MaCh3Exception(__FILE__ , __LINE__ );
278  }
279  if(BadXCounter == 5) MACH3LOG_WARN("There is more unutilised knots although I will stop spamming");
280  }
281 
282  MACH3LOG_WARN("Found in total {} BAD X", BadXCounter);
283  #ifdef Weight_On_SplineBySpline_Basis
284  // Make the array that holds all the returned weights from the GPU to pass to the CPU
285  cpu_weights_spline_var = new float[NSplines_valid]();
286  cpu_weights_tf1_var = new float[NTF1_valid]();
287  #else
288  //KS: This is tricky as this variable use both by CPU and GPU, however if use CUDA we use cudaMallocHost
289  #ifndef MaCh3_CUDA
290  cpu_total_weights = new float[NEvents]();
291  cpu_weights_spline_var = new float[NSplines_valid]();
292  cpu_weights_tf1_var = new float[NTF1_valid]();
293  #endif
294  #endif
295 
296  // Print some info; could probably make this to a separate function
299 
300  MoveToGPU();
301 }
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
@ kTF1_red
Uses TF1_red for interpolation.
@ kTSpline3_red
Uses TSpline3_red for interpolation.
std::vector< short int > cpu_nPoints_arr
CPU arrays to hold number of points.
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 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.
void getTF1Coeff(TF1_red *&spl, int &nPoints, float *&coeffs)
CW: Gets the polynomial coefficients for TF1.
Definition: SplineBase.cpp:106
CW: A reduced TF1 class only. Only saves parameters for each TF1 and how many parameters each paramet...
KS: A reduced ResponseFunction Generic function used for evaluating weight.
CW: Reduced TSpline3 class.
int int_t
Definition: Core.h:38

◆ PrepareSplineFile()

void SMonolith::PrepareSplineFile ( std::string  FileName)
overridevirtual

KS: Prepare spline file that can be used for fast loading.

Implements SplineBase.

Definition at line 608 of file SplineMonolith.cpp.

608  {
609 // *****************************************
610  M3::AddPath(FileName);
611 
612  auto SplineFile = std::make_unique<TFile>(FileName.c_str(), "recreate");
613  TTree *Settings = new TTree("Settings", "Settings");
614  TTree *Monolith_TF1 = new TTree("Monolith_TF1", "Monolith_TF1");
615  TTree *XKnots = new TTree("XKnots", "XKnots");
616  TTree *EventInfo = new TTree("EventInfo", "EventInfo");
617 
618  unsigned int NEvents_temp = NEvents;
619  short int nParams_temp = nParams;
620  int _max_knots_temp = _max_knots;
621  unsigned int nKnots_temp = nKnots;
622  unsigned int NSplines_valid_temp = NSplines_valid;
623  unsigned int nTF1Valid_temp = NTF1_valid;
624  unsigned int nTF1coeff_temp = nTF1coeff;
625 
626  Settings->Branch("NEvents", &NEvents_temp, "NEvents/i");
627  Settings->Branch("nParams", &nParams_temp, "nParams/S");
628  Settings->Branch("_max_knots", &_max_knots_temp, "_max_knots/I");
629  Settings->Branch("nKnots", &nKnots_temp, "nKnots/i");
630  Settings->Branch("NSplines_valid", &NSplines_valid_temp, "NSplines_valid/i");
631  Settings->Branch("NTF1_valid", &nTF1Valid_temp, "NTF1_valid/i");
632  Settings->Branch("nTF1coeff", &nTF1coeff_temp, "nTF1coeff/i");
633 
634  Settings->Fill();
635 
636  SplineFile->cd();
637  Settings->Write();
638 
639  TTree *SplineTree = new TTree("SplineTree", "SplineTree");
640  // Create a branch for the SplineMonoStruct object
641  SplineTree->Branch("SplineObject", &cpu_spline_handler);
642  SplineTree->Fill();
643  SplineTree->Write();
644  delete SplineTree;
645 
646  float coeff_tf1 = 0.;
647  Monolith_TF1->Branch("cpu_coeff_TF1_many", &coeff_tf1, "cpu_coeff_TF1_many/F");
648  for(unsigned int i = 0; i < nTF1coeff; i++)
649  {
650  coeff_tf1 = cpu_coeff_TF1_many[i];
651  Monolith_TF1->Fill();
652  }
653  SplineFile->cd();
654  Monolith_TF1->Write();
655 
656  unsigned int nParamPerEvent = 0;
657  unsigned int nParamPerEvent_tf1 = 0;
658 
659  EventInfo->Branch("cpu_nParamPerEvent", &nParamPerEvent, "cpu_nParamPerEvent/i");
660  EventInfo->Branch("cpu_nParamPerEvent_tf1", &nParamPerEvent_tf1, "cpu_nParamPerEvent_tf1/i");
661 
662  for(unsigned int i = 0; i < 2*NEvents; i++)
663  {
664  nParamPerEvent = cpu_nParamPerEvent[i];
665  nParamPerEvent_tf1 = cpu_nParamPerEvent_tf1[i];
666  EventInfo->Fill();
667  }
668  SplineFile->cd();
669  EventInfo->Write();
670 
671  PrepareFastSplineInfoDir(SplineFile);
672 
673  delete Settings;
674  delete Monolith_TF1;
675  delete XKnots;
676  delete EventInfo;
677  SplineFile->Close();
678 }
void PrepareFastSplineInfoDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Fast Spline Info within SplineFile.
Definition: SplineBase.cpp:131

◆ PrintInitialsiation()

void SMonolith::PrintInitialsiation ( )
inlineprivate

KS: Print info about how much knots etc has been initialised.

Definition at line 938 of file SplineMonolith.cpp.

938  {
939 //*********************************************************
940  unsigned int event_size_max = _max_knots * nParams;
941 
942  MACH3LOG_INFO("--- INITIALISED Spline Monolith ---");
943  MACH3LOG_INFO("{} events with {} splines", NEvents, NSplines_valid);
944  MACH3LOG_INFO("On average {:.2f} splines per event ({}/{})", float(NSplines_valid)/float(NEvents), NSplines_valid, NEvents);
945  MACH3LOG_INFO("Size of x array = {:.4f} MB", double(sizeof(float)*event_size_max)/1.E6);
946  MACH3LOG_INFO("Size of coefficient (y,b,c,d) array = {:.2f} MB", double(sizeof(float)*nKnots*_nCoeff_)/1.E6);
947  MACH3LOG_INFO("Size of parameter # array = {:.2f} MB", double(sizeof(short int)*NSplines_valid)/1.E6);
948 
949  MACH3LOG_INFO("On average {:.2f} TF1 per event ({}/{})", float(NTF1_valid)/float(NEvents), NTF1_valid, NEvents);
950  MACH3LOG_INFO("Size of TF1 coefficient (a,b,c,d,e) array = {:.2f} MB", double(sizeof(float)*NTF1_valid*_nTF1Coeff_)/1.E6);
951 }

◆ retPointer()

const float* SMonolith::retPointer ( const int  event) const
inline

KS: Get pointer to total weight to make fit faster wrooom!

Parameters
eventName event number in used MC
Returns
Pointer to the total weight

Definition at line 40 of file SplineMonolith.h.

40 {return &cpu_total_weights[event];}

◆ ScanMasterSpline()

void SMonolith::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 
)
inlineprivate

CW: Function to scan through the MasterSpline of TSpline3.

Parameters
MasterSplineVector of TSpline3_red pointers which we strip back
NEventsNumber of MC events
MaxPointsMaximal number of knots per splines
numParamsTotal number of parameters
numKnotsTotal number of knots, which is sum of individual knots per each spline
nTF1_coeffNumber of TF1 coefficients in all TF1 objects
SplineTypeWhether object is TSpline3 or TF1
NSplinesValidTotal number of valid (not null) TSpline3
nTF1ValidTotal number of valid (not null) TF1

Definition at line 377 of file SplineMonolith.cpp.

386  {
387 // *****************************************
388  // Need to extract: the total number of events
389  // number of parameters
390  // maximum number of knots
391  MaxPoints = 0;
392  nEvents = 0;
393  numParams = 0;
394  nSplines = 0;
395  numKnots = 0;
396  NSplinesValid = 0;
397  nTF1Valid = 0;
398  nTF1_coeff = 0;
399 
400  // Check the number of events
401  nEvents = int(MasterSpline.size());
402 
403  // Maximum number of splines one event can have (scan through and find this number)
404  int nMaxSplines_PerEvent = 0;
405 
406  //KS: We later check that each event has the same number of splines so this is fine
407  numParams = short(MasterSpline[0].size());
408  // Initialise
409  SplineInfoArray.resize(numParams);
410 
411  // Loop over each parameter
412  for(unsigned int EventCounter = 0; EventCounter < MasterSpline.size(); ++EventCounter) {
413  // Check that each event has each spline saved
414  if (numParams > 0) {
415  int TempSize = int(MasterSpline[EventCounter].size());
416  if (TempSize != numParams) {
417  MACH3LOG_ERROR("Found {} parameters for event {}", TempSize, EventCounter);
418  MACH3LOG_ERROR("but was expecting {} since that's what I found for the previous event", numParams);
419  MACH3LOG_ERROR("Somehow this event has a different number of spline parameters... Please study further!");
420  throw MaCh3Exception(__FILE__ , __LINE__ );
421  }
422  }
423  numParams = short(MasterSpline[EventCounter].size());
424 
425  int nSplines_SingleEvent = 0;
426  int nPoints = 0;
427  // Loop over each pointer
428  for(unsigned int ParamNumber = 0; ParamNumber < MasterSpline[EventCounter].size(); ++ParamNumber) {
429  if (MasterSpline[EventCounter][ParamNumber]) {
430  if(SplineType[ParamNumber] == kTSpline3_red)
431  {
432  TResponseFunction_red* TespFunc = MasterSpline[EventCounter][ParamNumber];
433  TSpline3_red* CurrSpline = dynamic_cast<TSpline3_red*>(TespFunc);
434  if(CurrSpline){
435  nPoints = CurrSpline->GetNp();
436  }
437 
438  if (nPoints > MaxPoints) {
439  MaxPoints = static_cast<short int>(nPoints);
440  }
441  numKnots += nPoints;
442  nSplines_SingleEvent++;
443 
444  // Fill the SplineInfoArray entries with information on each splinified parameter
445  if (SplineInfoArray[ParamNumber].xPts.size() == 0)
446  {
447  // Fill the number of points
448  SplineInfoArray[ParamNumber].nPts = CurrSpline->GetNp();
449 
450  // Fill the x points
451  SplineInfoArray[ParamNumber].xPts.resize(SplineInfoArray[ParamNumber].nPts);
452  for (M3::int_t k = 0; k < SplineInfoArray[ParamNumber].nPts; ++k)
453  {
454  M3::float_t xtemp = M3::float_t(-999.99);
455  M3::float_t ytemp = M3::float_t(-999.99);
456  CurrSpline->GetKnot(k, xtemp, ytemp);
457  SplineInfoArray[ParamNumber].xPts[k] = xtemp;
458  }
459  }
460  NSplinesValid++;
461  }
462  else if (SplineType[ParamNumber] == kTF1_red)
463  {
464  TResponseFunction_red* TespFunc = MasterSpline[EventCounter][ParamNumber];
465  TF1_red* CurrSpline = dynamic_cast<TF1_red*>(TespFunc);
466  nPoints = CurrSpline->GetSize();
467  nTF1_coeff += nPoints;
468  nTF1Valid++;
469  }
470  } else {
471  // If NULL we don't have this spline for the event, so move to next spline
472  continue;
473  }
474  }
475  if (nSplines_SingleEvent > nMaxSplines_PerEvent) nMaxSplines_PerEvent = nSplines_SingleEvent;
476  }
477  nSplines = nMaxSplines_PerEvent;
478 
479  int Counter = 0;
480  //KS: Sanity check that everything was set correctly
481  for (M3::int_t i = 0; i < numParams; ++i)
482  {
483  // KS: We don't find segment for TF1, so ignore this
484  if (SplineType[i] == kTF1_red) continue;
485 
486  const M3::int_t nPoints = SplineInfoArray[i].nPts;
487  const std::vector<M3::float_t>& xArray = SplineInfoArray[i].xPts;
488  if (nPoints == -999 || xArray.size() == 0) {
489  Counter++;
490  if(Counter < 5) {
491  MACH3LOG_WARN("SplineInfoArray[{}] isn't set yet", i);
492  }
493  continue;
494  //throw MaCh3Exception(__FILE__ , __LINE__ );
495  }
496  }
497  MACH3LOG_WARN("In total SplineInfoArray for {} hasn't been initialised", Counter);
498 }
std::vector< FastSplineInfo > SplineInfoArray
Definition: SplineBase.h:74
int GetSize()
Get the size.
void GetKnot(int i, M3::float_t &xtmp, M3::float_t &ytmp)

◆ setSplinePointers()

void SMonolith::setSplinePointers ( std::vector< const double * >  spline_ParsPointers)
inline

KS: Set pointers to spline params.

Parameters
spline_ParsPointersVector of pointers to spline params

Definition at line 44 of file SplineMonolith.h.

44  {
45  splineParsPointer = spline_ParsPointers;
46  for (M3::int_t i = 0; i < nParams; ++i) SplineInfoArray[i].splineParsPointer = spline_ParsPointers[i];
47  };
std::vector< const double * > splineParsPointer
This holds pointer to parameter position which we later copy paste it to GPU.

◆ SynchroniseMemTransfer()

void SMonolith::SynchroniseMemTransfer ( ) const
overridevirtual

KS: After calculations are done on GPU we copy memory to CPU. This operation is asynchronous meaning while memory is being copied some operations are being carried. Memory must be copied before actual reweight. This function make sure all has been copied.

Implements SplineBase.

Definition at line 955 of file SplineMonolith.cpp.

955  {
956 //*********************************************************
957  #ifdef MaCh3_CUDA
959  #endif
960 }
__host__ void SynchroniseSplines()
Make sure all Cuda threads finished execution.

Member Data Documentation

◆ _max_knots

short int SMonolith::_max_knots
private

Max knots for production.

Definition at line 113 of file SplineMonolith.h.

◆ cpu_coeff_TF1_many

std::vector<float> SMonolith::cpu_coeff_TF1_many
private

CPU arrays to hold TF1 coefficients.

Definition at line 150 of file SplineMonolith.h.

◆ cpu_nParamPerEvent

std::vector<unsigned int> SMonolith::cpu_nParamPerEvent
private

KS: CPU 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 138 of file SplineMonolith.h.

◆ cpu_nParamPerEvent_tf1

std::vector<unsigned int> SMonolith::cpu_nParamPerEvent_tf1
private

KS: CPU 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 141 of file SplineMonolith.h.

◆ cpu_nPoints_arr

std::vector<short int> SMonolith::cpu_nPoints_arr
private

CPU arrays to hold number of points.

Definition at line 153 of file SplineMonolith.h.

◆ cpu_paramNo_TF1_arr

std::vector<short int> SMonolith::cpu_paramNo_TF1_arr
private

CW: CPU array with the number of points per spline (not per spline point!)

Definition at line 156 of file SplineMonolith.h.

◆ cpu_spline_handler

SplineMonoStruct* SMonolith::cpu_spline_handler
private

KS: Store info about Spline monolith, this allow to obtain better step time. As all necessary information for spline weight calculation are here meaning better cache hits.

Definition at line 144 of file SplineMonolith.h.

◆ cpu_total_weights

float* SMonolith::cpu_total_weights

KS: This holds the total CPU weights that gets read in samplePDFND.

Definition at line 52 of file SplineMonolith.h.

◆ cpu_weights

float* SMonolith::cpu_weights

The returned gpu weights, read by the GPU.

Definition at line 50 of file SplineMonolith.h.

◆ cpu_weights_spline_var

float* SMonolith::cpu_weights_spline_var
private

CPU arrays to hold weight for each spline.

Definition at line 133 of file SplineMonolith.h.

◆ cpu_weights_tf1_var

float* SMonolith::cpu_weights_tf1_var
private

CPU arrays to hold weight for each TF1.

Definition at line 135 of file SplineMonolith.h.

◆ FastSplineName

std::string SMonolith::FastSplineName
private

Name of Fast Spline to which will be saved.

Definition at line 162 of file SplineMonolith.h.

◆ gpu_spline_handler

SMonolithGPU* SMonolith::gpu_spline_handler
private

KS: Store info about Spline monolith, this allow to obtain better step time. As all necessary information for spline weight calculation are here meaning better cache hits.

Definition at line 147 of file SplineMonolith.h.

◆ index_spline_cpu

std::vector<int> SMonolith::index_spline_cpu
private

holds the index for good splines; don't do unsigned since starts with negative value!

Definition at line 115 of file SplineMonolith.h.

◆ index_TF1_cpu

std::vector<int> SMonolith::index_TF1_cpu
private

holds the index for good TF1; don't do unsigned since starts with negative value!

Definition at line 117 of file SplineMonolith.h.

◆ NEvents

unsigned int SMonolith::NEvents
private

Number of events.

Definition at line 111 of file SplineMonolith.h.

◆ nKnots

unsigned int SMonolith::nKnots
private

Sum of all knots over all splines.

Definition at line 128 of file SplineMonolith.h.

◆ NSplines_total_large

unsigned int SMonolith::NSplines_total_large
private

Number of total splines if each event had every parameter's spline.

Definition at line 125 of file SplineMonolith.h.

◆ NSplines_valid

unsigned int SMonolith::NSplines_valid
private

Number of valid splines.

Definition at line 120 of file SplineMonolith.h.

◆ NTF1_valid

unsigned int SMonolith::NTF1_valid
private

Number of valid TF1.

Definition at line 122 of file SplineMonolith.h.

◆ nTF1coeff

unsigned int SMonolith::nTF1coeff
private

Sum of all coefficients over all TF1.

Definition at line 130 of file SplineMonolith.h.

◆ SaveSplineFile

bool SMonolith::SaveSplineFile
private

Flag telling whether we are saving spline monolith into handy root file.

Definition at line 159 of file SplineMonolith.h.

◆ splineParsPointer

std::vector< const double* > SMonolith::splineParsPointer
private

This holds pointer to parameter position which we later copy paste it to GPU.

Definition at line 108 of file SplineMonolith.h.


The documentation for this class was generated from the following files: