MaCh3  2.5.0
Reference Guide
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
SplineBase Class Referenceabstract

Base class for calculating weight from spline. More...

#include <Splines/SplineBase.h>

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

Public Member Functions

 SplineBase ()
 Constructor. More...
 
virtual ~SplineBase ()
 Destructor. More...
 
virtual void Evaluate ()=0
 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...
 
virtual std::string GetName () const
 Get class name. More...
 
short int GetNParams () const
 Get number of spline parameters. More...
 
virtual void PrepareSplineFile (std::string FileName)=0
 KS: Prepare spline file that can be used for fast loading. More...
 
virtual void LoadSplineFile (std::string FileName)=0
 KS: Load preprocessed spline file. More...
 
virtual void SynchroniseMemTransfer () const =0
 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...
 

Protected Member Functions

void FindSplineSegment ()
 CW:Code used in step by step reweighting, Find Spline Segment for each param. More...
 
virtual void CalcSplineWeights ()=0
 CPU based code which eval weight for each spline. 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) const
 CW: Gets the polynomial coefficients for TF1. More...
 

Protected Attributes

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...
 

Private Member Functions

void CheckSegmentValidity (const M3::int_t i, const M3::int_t segment, const M3::float_t xvar, const std::vector< M3::float_t > &xPts) const
 Validates that the spline segment is correct for the given variation. More...
 

Detailed Description

Base class for calculating weight from spline.

Author
Dan Barrow
Kamil Skwarczynski

Definition at line 27 of file SplineBase.h.

Constructor & Destructor Documentation

◆ SplineBase()

SplineBase::SplineBase ( )

Constructor.

Definition at line 7 of file SplineBase.cpp.

7  {
8 // *****************************************
9  nParams = 0;
10  SplineSegments = nullptr;
11  ParamValues = nullptr;
12 }
short int nParams
Number of parameters that have splines.
Definition: SplineBase.h:81
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

◆ ~SplineBase()

SplineBase::~SplineBase ( )
virtual

Destructor.

Definition at line 16 of file SplineBase.cpp.

16  {
17 // *****************************************
18 }

Member Function Documentation

◆ CalcSplineWeights()

virtual void SplineBase::CalcSplineWeights ( )
protectedpure virtual

CPU based code which eval weight for each spline.

Implemented in PySplineBase, UnbinnedSplineHandler, and BinnedSplineHandler.

◆ CheckSegmentValidity()

void SplineBase::CheckSegmentValidity ( const M3::int_t  i,
const M3::int_t  segment,
const M3::float_t  xvar,
const std::vector< M3::float_t > &  xPts 
) const
private

Validates that the spline segment is correct for the given variation.

Parameters
iIndex of the spline in SplineInfoArray.
segmentThe segment index to validate.
xvarThe variation value being checked.
xPtsThe array of x-values for the spline.

Definition at line 21 of file SplineBase.cpp.

24  {
25 // *************************
26  if (xPts[segment] > xvar && segment != 0) {
27  MACH3LOG_ERROR("Found a segment which is _ABOVE_ the variation!");
28  MACH3LOG_ERROR("IT SHOULD ALWAYS BE BELOW! (except when segment 0)");
29  MACH3LOG_ERROR("Spline: {}", i);
30  MACH3LOG_ERROR("Found segment = {}", segment);
31  MACH3LOG_ERROR("Doing variation = {}", xvar);
32  MACH3LOG_ERROR("x in spline = {}", xPts[segment]);
33  for (size_t j = 0; j < xPts.size(); ++j) {
34  MACH3LOG_ERROR(" {} = {}", j, xPts[j]);
35  }
36  throw MaCh3Exception(__FILE__, __LINE__);
37  }
38 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
Custom exception class used throughout MaCh3.

◆ Evaluate()

virtual void SplineBase::Evaluate ( )
pure virtual

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.

Implemented in PySplineBase, UnbinnedSplineHandler, and BinnedSplineHandler.

◆ FindSplineSegment()

void SplineBase::FindSplineSegment ( )
protected

CW:Code used in step by step reweighting, Find Spline Segment for each param.

Definition at line 44 of file SplineBase.cpp.

44  {
45 // *************************
46  // Loop over the splines
47  //KS: Tried multithreading here with 48 splines and it is faster with one thread, maybe in future multithreading will be worth revisiting
48  for (M3::int_t i = 0; i < nParams; ++i)
49  {
50  const M3::int_t nPoints = SplineInfoArray[i].nPts;
51  const std::vector<M3::float_t>& xArray = SplineInfoArray[i].xPts;
52 
53  // Get the variation for this reconfigure for the ith parameter
54  const float xvar = float(*SplineInfoArray[i].splineParsPointer);
55  ParamValues[i] = xvar;
56 
57  // EM: if we have a parameter that has no response for any event (i.e. all splines have just one knot), then skip it and avoid a seg fault here
58  // In principle, such parameters shouldn't really be included in the first place, but with new det syst splines this
59  // could happen if say you were just running on one FHC run, then all RHC parameters would be flat and the code below would break.
60  if(xArray.size() == 0) continue;
61 
62  // The segment we're interested in (klow in ROOT code)
63  M3::int_t segment = 0;
64  M3::int_t kHigh = nPoints-1;
65  //KS: We expect new segment is very close to previous
66  const M3::int_t PreviousSegment = SplineInfoArray[i].CurrSegment;
67 
68  // If the variation is below the lowest saved spline point
69  if (xvar <= xArray[0]) {
70  segment = 0;
71  // If the variation is above the highest saved spline point
72  } else if (xvar >= xArray[nPoints-1]) {
73  //CW: Yes, the -2 is indeed correct, see TSpline.cxx:814 and //see: https://savannah.cern.ch/bugs/?71651
74  segment = kHigh;
75  //KS: It is quite probable the new segment is same as in previous step so try to avoid binary search, first we have to check if it is in bounds to avoid seg fault
76  } else if( xArray[PreviousSegment+1] > xvar && xvar >= xArray[PreviousSegment] ) {
77  segment = PreviousSegment;
78  // If the variation is between the maximum and minimum, perform a binary search
79  } else {
80  // The top point we've got
81  M3::int_t kHalf = 0;
82  // While there is still a difference in the points (we haven't yet found the segment)
83  // This is a binary search, incrementing segment and decrementing kHalf until we've found the segment
84  while (kHigh - segment > 1) {
85  // Increment the half-step
86  kHalf = M3::int_t((segment + kHigh)/2);
87  // If our variation is above the kHalf, set the segment to kHalf
88  if (xvar > xArray[kHalf]) {
89  segment = kHalf;
90  // Else move kHigh down
91  } else {
92  kHigh = kHalf;
93  }
94  } // End the while: we've now done our binary search
95  } // End the else: we've now found our point
96 
97  if (segment >= nPoints-1 && nPoints > 1) segment = nPoints-2;
98 
99  // CW: This way we avoid doing 1.2M+ binary searches on the GPU
100  // and literally just multiply lots of numbers together on the GPU without any algorithm
101  // Update the values and which segment it belongs to
102  SplineInfoArray[i].CurrSegment = segment;
103  SplineSegments[i] = short(SplineInfoArray[i].CurrSegment);
104 
105  #ifdef MACH3_DEBUG
106  CheckSegmentValidity(i, segment, xvar, SplineInfoArray[i].xPts);
107  #endif
108  } //end loop over params
109 }
void CheckSegmentValidity(const M3::int_t i, const M3::int_t segment, const M3::float_t xvar, const std::vector< M3::float_t > &xPts) const
Validates that the spline segment is correct for the given variation.
Definition: SplineBase.cpp:21
std::vector< FastSplineInfo > SplineInfoArray
Definition: SplineBase.h:74
int int_t
Definition: Core.h:38

◆ GetName()

virtual std::string SplineBase::GetName ( ) const
inlinevirtual

Get class name.

Reimplemented in UnbinnedSplineHandler, and PySplineBase.

Definition at line 38 of file SplineBase.h.

38 {return "SplineBase";};

◆ GetNParams()

short int SplineBase::GetNParams ( ) const
inline

Get number of spline parameters.

Definition at line 41 of file SplineBase.h.

41 {return nParams;};

◆ GetTF1Coeff()

void SplineBase::GetTF1Coeff ( TF1_red *&  spl,
int &  nPoints,
float *&  coeffs 
) const
protected

CW: Gets the polynomial coefficients for TF1.

Parameters
splpointer to TF1_red that will be checked
nPointsnumber of knots
coeffsArray holding coefficients for each knot

Definition at line 115 of file SplineBase.cpp.

115  {
116 // *****************************************
117  // Initialise all arrays to 1.0
118  for (int i = 0; i < _nTF1Coeff_; ++i) {
119  coeffs[i] = 0.0;
120  }
121 
122  // Get number of points in spline
123  nPoints = spl->GetSize();
124 
125  if(nPoints > _nTF1Coeff_)
126  {
127  MACH3LOG_ERROR("Too big number of points for TF1");
128  throw MaCh3Exception(__FILE__ , __LINE__ );
129  }
130  // TSpline3 can only take doubles, not floats
131  // But our GPU is slow with doubles, so need to cast to float
132  for (int i = 0; i < nPoints; i++) {
133  coeffs[i] = float(spl->GetParameter(M3::int_t(i)));
134  }
135  // The structure is now coeffs = {a,b,c,d,e}
136 }
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:20
double GetParameter(M3::int_t Parameter)
Get a parameter value.
int GetSize()
Get the size.

◆ LoadFastSplineInfoDir()

void SplineBase::LoadFastSplineInfoDir ( std::unique_ptr< TFile > &  SplineFile)
protected

KS: Load preprocessed FastSplineInfo.

Parameters
FileFile from which we load new tree

Definition at line 167 of file SplineBase.cpp.

167  {
168 // *****************************************
169  TTree *FastSplineInfoTree = SplineFile->Get<TTree>("FastSplineInfoTree");
170  M3::int_t nPoints = 0;
171  float xtemp[20];
172  FastSplineInfoTree->SetBranchAddress("nPts", &nPoints);
173  FastSplineInfoTree->SetBranchAddress("xPts", &xtemp);
174 
175  if(nParams == 0){
176  MACH3LOG_WARN("Calling {} with {} number of params", __func__, nParams);
177  }
178  SplineInfoArray.resize(nParams);
179  for (M3::int_t i = 0; i < nParams; ++i) {
180  FastSplineInfoTree->GetEntry(i);
181 
182  // Fill the number of points
183  SplineInfoArray[i].nPts = nPoints;
184  if(nPoints == -999) continue;
185  SplineInfoArray[i].xPts.resize(SplineInfoArray[i].nPts);
186  for (M3::int_t k = 0; k < SplineInfoArray[i].nPts; ++k)
187  {
188  SplineInfoArray[i].xPts[k] = xtemp[k];
189  }
190  }
191 }
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36

◆ LoadSplineFile()

virtual void SplineBase::LoadSplineFile ( std::string  FileName)
pure virtual

KS: Load preprocessed spline file.

Parameters
FileNamePath to ROOT file with predefined reduced Spline Monolith

Implemented in UnbinnedSplineHandler, and BinnedSplineHandler.

◆ PrepareFastSplineInfoDir()

void SplineBase::PrepareFastSplineInfoDir ( std::unique_ptr< TFile > &  SplineFile) const
protected

KS: Prepare Fast Spline Info within SplineFile.

Parameters
FileFile to which we add new tree

Definition at line 139 of file SplineBase.cpp.

139  {
140 // *****************************************
141  TTree *FastSplineInfoTree = new TTree("FastSplineInfoTree", "FastSplineInfoTree");
142 
143  M3::int_t nPoints = 0;
144  float xtemp[20];
145  FastSplineInfoTree->Branch("nPts", &nPoints, "nPts/I");
146  FastSplineInfoTree->Branch("xPts", xtemp, "xPts[nPts]/F");
147 
148  for (M3::int_t i = 0; i < nParams; ++i)
149  {
150  nPoints = SplineInfoArray[i].nPts;
151 
152  for (M3::int_t k = 0; k < SplineInfoArray[i].nPts; ++k)
153  {
154  xtemp[k] = float(SplineInfoArray[i].xPts[k]);
155  }
156  FastSplineInfoTree->Fill();
157  }
158 
159  SplineFile->cd();
160  FastSplineInfoTree->Write();
161 
162  delete FastSplineInfoTree;
163 }

◆ PrepareSplineFile()

virtual void SplineBase::PrepareSplineFile ( std::string  FileName)
pure virtual

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

Implemented in UnbinnedSplineHandler, and BinnedSplineHandler.

◆ SynchroniseMemTransfer()

virtual void SplineBase::SynchroniseMemTransfer ( ) const
pure virtual

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.

Implemented in UnbinnedSplineHandler, and BinnedSplineHandler.

Member Data Documentation

◆ nParams

short int SplineBase::nParams
protected

Number of parameters that have splines.

Definition at line 81 of file SplineBase.h.

◆ ParamValues

float* SplineBase::ParamValues
protected

Store parameter values they are not in FastSplineInfo as in case of GPU we need to copy paste it to GPU.

Definition at line 79 of file SplineBase.h.

◆ SplineInfoArray

std::vector<FastSplineInfo> SplineBase::SplineInfoArray
protected

Array of FastSplineInfo structs: keeps information on each xsec spline for fast evaluation Method identical to TSpline3::Eval(double) but faster because less operations

Definition at line 74 of file SplineBase.h.

◆ SplineSegments

short int* SplineBase::SplineSegments
protected

Store currently found segment they are not in FastSplineInfo as in case of GPU we need to copy paste it to GPU

Warning
this is being used sometimes by GPU, therefore must be raw pointer!

Definition at line 77 of file SplineBase.h.


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