Base class for calculating weight from spline.
More...
#include <Splines/SplineBase.h>
|
| | 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...
|
| |
Base class for calculating weight from spline.
- Author
- Dan Barrow
-
Kamil Skwarczynski
Definition at line 27 of file SplineBase.h.
◆ SplineBase()
| SplineBase::SplineBase |
( |
| ) |
|
Constructor.
Definition at line 7 of file SplineBase.cpp.
short int nParams
Number of parameters that have splines.
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 G...
◆ ~SplineBase()
| SplineBase::~SplineBase |
( |
| ) |
|
|
virtual |
◆ CalcSplineWeights()
| virtual void SplineBase::CalcSplineWeights |
( |
| ) |
|
|
protectedpure virtual |
◆ CheckSegmentValidity()
Validates that the spline segment is correct for the given variation.
- Parameters
-
| i | Index of the spline in SplineInfoArray. |
| segment | The segment index to validate. |
| xvar | The variation value being checked. |
| xPts | The array of x-values for the spline. |
Definition at line 21 of file SplineBase.cpp.
26 if (xPts[segment] > xvar && segment != 0) {
28 MACH3LOG_ERROR(
"IT SHOULD ALWAYS BE BELOW! (except when segment 0)");
33 for (
size_t j = 0; j < xPts.size(); ++j) {
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.
60 if(xArray.size() == 0)
continue;
69 if (xvar <= xArray[0]) {
72 }
else if (xvar >= xArray[nPoints-1]) {
76 }
else if( xArray[PreviousSegment+1] > xvar && xvar >= xArray[PreviousSegment] ) {
77 segment = PreviousSegment;
84 while (kHigh - segment > 1) {
88 if (xvar > xArray[kHalf]) {
97 if (segment >= nPoints-1 && nPoints > 1) segment = nPoints-2;
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.
std::vector< FastSplineInfo > SplineInfoArray
◆ GetName()
| virtual std::string SplineBase::GetName |
( |
| ) |
const |
|
inlinevirtual |
◆ GetNParams()
| short int SplineBase::GetNParams |
( |
| ) |
const |
|
inline |
Get number of spline parameters.
Definition at line 41 of file SplineBase.h.
◆ GetTF1Coeff()
| void SplineBase::GetTF1Coeff |
( |
TF1_red *& |
spl, |
|
|
int & |
nPoints, |
|
|
float *& |
coeffs |
|
) |
| const |
|
protected |
CW: Gets the polynomial coefficients for TF1.
- Parameters
-
| spl | pointer to TF1_red that will be checked |
| nPoints | number of knots |
| coeffs | Array holding coefficients for each knot |
Definition at line 115 of file SplineBase.cpp.
132 for (
int i = 0; i < nPoints; i++) {
constexpr int _nTF1Coeff_
KS: For TF1 we store at most 5 coefficients, we could make it more flexible but for now define it her...
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
-
| File | File from which we load new tree |
Definition at line 167 of file SplineBase.cpp.
169 TTree *FastSplineInfoTree = SplineFile->Get<TTree>(
"FastSplineInfoTree");
172 FastSplineInfoTree->SetBranchAddress(
"nPts", &nPoints);
173 FastSplineInfoTree->SetBranchAddress(
"xPts", &xtemp);
180 FastSplineInfoTree->GetEntry(i);
184 if(nPoints == -999)
continue;
◆ LoadSplineFile()
| virtual void SplineBase::LoadSplineFile |
( |
std::string |
FileName | ) |
|
|
pure virtual |
◆ PrepareFastSplineInfoDir()
| void SplineBase::PrepareFastSplineInfoDir |
( |
std::unique_ptr< TFile > & |
SplineFile | ) |
const |
|
protected |
KS: Prepare Fast Spline Info within SplineFile.
- Parameters
-
| File | File to which we add new tree |
Definition at line 139 of file SplineBase.cpp.
141 TTree *FastSplineInfoTree =
new TTree(
"FastSplineInfoTree",
"FastSplineInfoTree");
145 FastSplineInfoTree->Branch(
"nPts", &nPoints,
"nPts/I");
146 FastSplineInfoTree->Branch(
"xPts", xtemp,
"xPts[nPts]/F");
156 FastSplineInfoTree->Fill();
160 FastSplineInfoTree->Write();
162 delete FastSplineInfoTree;
◆ PrepareSplineFile()
| virtual void SplineBase::PrepareSplineFile |
( |
std::string |
FileName | ) |
|
|
pure virtual |
◆ 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.
◆ 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
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: