![]() |
MaCh3 2.2.1
Reference Guide
|
Base class for calculating weight from spline. More...
#include <Splines/SplineBase.h>
Public Member Functions | |
SplineBase () | |
Constructor. | |
virtual | ~SplineBase () |
Destructor. | |
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. | |
virtual std::string | GetName () const |
Get class name. | |
short int | GetNParams () const |
Get number of spline parameters. | |
Protected Member Functions | |
void | FindSplineSegment () |
CW:Code used in step by step reweighting, Find Spline Segment for each param. | |
virtual void | CalcSplineWeights ()=0 |
CPU based code which eval weight for each spline. | |
virtual void | ModifyWeights ()=0 |
Calc total event weight. | |
void | getTF1Coeff (TF1_red *&spl, int &nPoints, float *&coeffs) |
CW: Gets the polynomial coefficients for TF1. | |
Protected Attributes | |
std::vector< FastSplineInfo > | SplineInfoArray |
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. | |
short int | nParams |
Number of parameters that have splines. | |
Base class for calculating weight from spline.
Definition at line 25 of file SplineBase.h.
SplineBase::SplineBase | ( | ) |
Constructor.
Definition at line 7 of file SplineBase.cpp.
|
virtual |
Destructor.
Definition at line 16 of file SplineBase.cpp.
|
protectedpure virtual |
CPU based code which eval weight for each spline.
Implemented in PySplineBase, BinnedSplineHandler, and SMonolith.
|
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 BinnedSplineHandler, PySplineBase, and SMonolith.
|
protected |
CW:Code used in step by step reweighting, Find Spline Segment for each param.
Definition at line 24 of file SplineBase.cpp.
|
inlinevirtual |
Get class name.
Reimplemented in SMonolith, and PySplineBase.
Definition at line 36 of file SplineBase.h.
|
inline |
|
protected |
CW: Gets the polynomial coefficients for TF1.
spl | pointer to TF1_red that will be checked |
nPoints | number of knots |
coeffs | Array holding coefficients for each knot |
Definition at line 106 of file SplineBase.cpp.
|
protectedpure virtual |
Calc total event weight.
Implemented in PySplineBase, BinnedSplineHandler, and SMonolith.
|
protected |
Number of parameters that have splines.
Definition at line 64 of file SplineBase.h.
|
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 62 of file SplineBase.h.
|
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 57 of file SplineBase.h.
|
protected |
Store currently found segment they are not in FastSplineInfo as in case of GPU we need to copy paste it to GPU
Definition at line 60 of file SplineBase.h.