9 #pragma GCC diagnostic push
10 #pragma GCC diagnostic ignored "-Wuseless-cast"
11 #pragma GCC diagnostic ignored "-Wfloat-conversion"
12 #pragma GCC diagnostic ignored "-Wnull-dereference"
37 std::vector<M3::float_t>
xPts;
51 if(xsecgraph ==
nullptr){
52 MACH3LOG_ERROR(
"hmmm looks like you're trying to apply capping for spline parameter {} but it hasn't been set in xsecgraph yet", XsecCov->
GetParFancyName(splineParsIndex));
61 for(
int knotId = 0; knotId < xsecgraph->GetN(); knotId++){
65 if( xsecgraph->GetPoint(knotId, x, y) == -1) {
73 xsecgraph->SetPoint(knotId, x, y);
78 for(
int knotId = 0; knotId < xsecgraph->GetN(); knotId++){
82 if( xsecgraph->GetPoint(knotId, x, y) == -1) {
86 if(std::abs(y - 1.0) > 1e-5)
isFlat =
false;
91 xsecgraph->SetPoint(0, 0.0, 1.0);
100 if(Spline ==
nullptr) {
101 MACH3LOG_ERROR(
"hmmm looks like you're trying to apply capping for spline parameter {} but it hasn't been set in Spline yet", XsecCov->
GetParFancyName(splineParsIndex));
105 std::string oldName = Spline->GetName();
109 const int NValues = Spline->GetNp();
110 std::vector<double> XVals(NValues);
111 std::vector<double> YVals(NValues);
112 for(
int knotId = 0; knotId < NValues; knotId++){
115 Spline->GetKnot(knotId, x, y);
125 Spline =
new TSpline3(oldName.c_str(), XVals.data(), YVals.data(), NValues);
139 virtual double Eval(
const double var)=0;
168 for (
int i = 0; i <
length; ++i) {
182 if (
Par !=
nullptr)
delete[]
Par;
184 for (
int i = 0; i <
length; ++i) {
192 inline double Eval(
const double var)
override {
216 Par[Parameter] = Value;
222 MACH3LOG_ERROR(
"You requested parameter number {} but length is {} parameters", Parameter,
length);
226 return Par[Parameter];
240 for (
int i = 0; i <
length; i++) {
246 inline TF1*
ConstructTF1(
const std::string&
function,
const int xmin,
const int xmax) {
247 TF1 *func =
new TF1(
"TF1",
function.c_str(), xmin, xmax);
248 for(
int i = 0; i <
length; ++i) {
249 func->SetParameter(i,
Par[i]);
281 SetFunc(spline, InterPolation);
293 for(
int j = 0; j < N; ++j){
301 if((
Par[j][0] == -999) | (
Par[j][1] ==-999) | (
Par[j][2] ==-999) | (
XPos[j] ==-999) | (
YResp[j] ==-999)){
302 MACH3LOG_ERROR(
"******************* Bad parameter values when constructing TSpline3_red *********************");
303 MACH3LOG_ERROR(
"Passed values (i, x, y, b, c, d): {}, {}, {}, {}, {}, {}", j, X[j], Y[j], P[j][0], P[j][1], P[j][2]);
304 MACH3LOG_ERROR(
"Set values (i, x, y, b, c, d): {}, {}, {}, {}, {}, {}", j,
XPos[j],
YResp[j],
Par[j][0],
Par[j][1],
Par[j][2]);
305 MACH3LOG_ERROR(
"*********************************************************************************************");
313 if (
Par !=
nullptr) {
314 for (
int i = 0; i <
nPoints; ++i) {
333 for (
int i = 0; i <
nPoints; ++i) {
336 double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
337 spline->GetCoeff(i, x, y, b, c, d);
353 for (
int k = 0; k <
nPoints; ++k) {
356 Double_t x1, y1, b1, c1, d1, x2, y2, b2, c2, d2 = 0;
357 spline->GetCoeff(k, x1, y1, b1, c1, d1);
363 spline->GetCoeff(k + 1, x2, y2, b2, c2, d2);
364 tempb = (y2-y1)/(x2-x1);
375 else if(InterPolation ==
kAkima)
378 for (
int i = 0; i <
nPoints; ++i) {
382 double x = -999.99, y = -999.99;
383 spline->GetKnot(i, x, y);
392 for (
int i = -2; i <=
nPoints; ++i) {
413 for(
int i =2; i<=
nPoints+2; i++){
414 if (std::abs(mvals[i+1] - mvals[i]) + std::abs(mvals[i-1] - mvals[i-2]) != 0.0){
415 svals[i-2] = (std::abs(mvals[i+1] - mvals[i]) * mvals[i-1] + std::abs(mvals[i-1] - mvals[i-2]) *mvals[i]) / (std::abs(mvals[i+1] - mvals[i]) + std::abs(mvals[i-1] - mvals[i-2]));
417 else{svals[i-2] = mvals[i];}
421 for(
int i = 0; i <
nPoints; i++){
435 for(
int i = 0; i <
nPoints-1; i++){
436 double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
437 spline->GetCoeff(i, x, y, b, c, d);
439 if((c == 0.0 && d == 0.0)){
458 for (
int i = 0; i <
nPoints; ++i) {
462 double x = -999.99, y = -999.99;
463 spline->GetKnot(i, x, y);
485 for (
int i = 0; i <
nPoints-1; ++i) {
490 Tangents[0] = Secants[0];
497 for (
int i = 1; i <
nPoints-1; ++i) {
498 if ((Secants[i-1] >= 0.0 && Secants[i] >= 0.0) | (Secants[i-1] < 0.0 && Secants[i] < 0.0)){
499 Tangents[i] =
M3::float_t((Secants[i-1] + Secants[i]) /2.0);
504 for (
int i = 0; i <
nPoints-1; ++i) {
505 if (Secants[i] == 0.0){
511 alpha = Tangents[i] / Secants[i];
512 beta = Tangents[i+1] / Secants[i];
521 if (alpha * alpha + beta * beta >9.0){
523 Tangents[i] = tau * alpha * Secants[i];
524 Tangents[i+1] = tau * beta * Secants[i];
530 for(
int i = 0; i <
nPoints-1; i++){
534 b = Tangents[i] * dx;
539 Par[i][1] = c / (dx * dx);
540 Par[i][2] = d / (dx * dx * dx);
542 if((
Par[i][0] == -999) | (
Par[i][1] == -999) | (
Par[i][2] ==-999) | (
Par[i][0] == -999.999) | (
Par[i][1] == -999.999) | (
Par[i][2] ==-999.999)){
543 MACH3LOG_INFO(
"Bad spline parameters for segment {}: (b, c, d) = {}, {}, {}. This will cause problems with GPU.",
559 for(
int i = 0; i <
nPoints-1; i++){
560 double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
561 spline->GetCoeff(i, x, y, b, c, d);
563 if((c == 0.0 && d == 0.0)){
577 for (
int i = 0; i <
nPoints; ++i)
580 double x = -999.99, y = -999.99;
581 spline->GetKnot(i, x, y);
607 for (
int i = 1; i <
nPoints - 1; ++i)
612 M3::float_t term1 = (1.0 - T) * (1.0 + B) * (1.0 + C) * 0.5 * d0;
613 M3::float_t term2 = (1.0 - T) * (1.0 - B) * (1.0 - C) * 0.5 * d1;
615 Tangents[i] = term1 + term2;
623 for (
int i = 0; i <
nPoints - 1; ++i)
631 Par[i][1] = (3*dy/(dx*dx)) - (2*m0 + m1)/dx;
632 Par[i][2] = (m0 + m1 - 2*dy/dx) / (dx*dx);
635 i, dx, dy,
Par[i][0],
Par[i][1],
Par[i][2]);
646 MACH3LOG_ERROR(
"Unsupported interpolation type {}",
static_cast<int>(InterPolation));
657 for (
int i = 0; i <
nPoints; ++i) {
658 if (
Par[i] !=
nullptr) {
690 while (kHigh - segment > 1) {
692 kHalf = (segment + kHigh)/2;
694 if (x >
XPos[kHalf]) {
707 inline double Eval(
double var)
override {
709 int segment =
FindX(var);
715 double weight = y+dx*(b+dx*(c+d*dx));
739 #ifdef _LOW_MEMORY_STRUCTS_
740 std::vector<Double_t> xPosDoubles(
nPoints);
741 std::vector<Double_t> yPosDoubles(
nPoints);
742 for (Int_t i = 0; i <
nPoints; ++i) {
743 xPosDoubles[i] =
static_cast<Double_t
>(
XPos[i]);
744 yPosDoubles[i] =
static_cast<Double_t
>(
YResp[i]);
746 TSpline3 *spline =
new TSpline3(
"Spline", xPosDoubles.data(), yPosDoubles.data(),
static_cast<int>(
nPoints));
750 for (Int_t i = 0; i <
nPoints; ++i) {
751 spline->SetPointCoeff(i,
Par[i][0],
Par[i][1],
Par[i][2]);
761 for (
int i = 0; i <
nPoints; ++i) {
783 int Np = spl->
GetNp();
788 for(
int i = 0; i < Np; i++) {
800 inline std::vector<std::vector<TSpline3_red*> >
ReduceTSpline3(std::vector<std::vector<TSpline3*> > &MasterSpline) {
802 std::vector<std::vector<TSpline3*> >::iterator OuterIt;
803 std::vector<TSpline3*>::iterator InnerIt;
806 std::vector<std::vector<TSpline3_red*> > ReducedVector;
807 ReducedVector.reserve(MasterSpline.size());
810 int OuterCounter = 0;
811 for (OuterIt = MasterSpline.begin(); OuterIt != MasterSpline.end(); ++OuterIt, ++OuterCounter) {
813 std::vector<TSpline3_red*> TempVector;
814 TempVector.reserve(OuterIt->size());
815 int InnerCounter = 0;
817 for (InnerIt = OuterIt->begin(); InnerIt != OuterIt->end(); ++InnerIt, ++InnerCounter) {
819 TSpline3 *spline = (*InnerIt);
822 if (spline != NULL) {
827 TempVector.push_back(red);
829 ReducedVector.push_back(TempVector);
832 return ReducedVector;
838 inline std::vector<std::vector<TF1_red*> >
ReduceTF1(std::vector<std::vector<TF1*> > &MasterSpline) {
840 std::vector<std::vector<TF1*> >::iterator OuterIt;
841 std::vector<TF1*>::iterator InnerIt;
844 std::vector<std::vector<TF1_red*> > ReducedVector;
845 ReducedVector.reserve(MasterSpline.size());
848 int OuterCounter = 0;
849 for (OuterIt = MasterSpline.begin(); OuterIt != MasterSpline.end(); ++OuterIt, ++OuterCounter) {
851 std::vector<TF1_red*> TempVector;
852 TempVector.reserve(OuterIt->size());
853 int InnerCounter = 0;
855 for (InnerIt = OuterIt->begin(); InnerIt != OuterIt->end(); ++InnerIt, ++InnerCounter) {
857 TF1* spline = (*InnerIt);
860 if (spline != NULL) {
865 TempVector.push_back(red);
867 ReducedVector.push_back(TempVector);
870 return ReducedVector;
882 const std::string& Title) {
886 if (graph && graph->GetN() > 1)
891 TSpline3* spline =
nullptr;
895 spline =
new TSpline3(Title.c_str(), graph);
896 spline->SetNameTitle(Title.c_str(), Title.c_str());
899 spline_red =
new TSpline3_red(spline, SplineInterpolationType);
901 RespFunc = spline_red;
903 else if(SplineRespFuncType ==
kTF1_red)
906 TF1 *Fitter =
nullptr;
909 if(graph->GetN() != 2) {
910 MACH3LOG_ERROR(
"Trying to make TF1 from more than 2 knots. Knots = {}", graph->GetN());
916 Fitter =
new TF1(Title.c_str(),
"([1]+[0]*x)", graph->GetX()[0], graph->GetX()[graph->GetN()-1]);
925 graph->Fit(Fitter,
"Q0");
944 #pragma GCC diagnostic pop
Definitions of generic parameter structs and utility templates for MaCh3.
SplineInterpolation
Make an enum of the spline interpolation type.
@ kTSpline3
Default TSpline3 interpolation.
@ kMonotonic
EM: DOES NOT make the entire spline monotonic, only the segments.
@ kKochanekBartels
KS: Kochanek-Bartels spline: allows local control of tension, continuity, and bias.
@ kLinear
Linear interpolation between knots.
@ kLinearFunc
Liner interpolation using TF1 not spline.
@ kAkima
EM: Akima spline iis allowed to be discontinuous in 2nd derivative and coefficients in any segment.
RespFuncType
Make an enum of the spline interpolation type.
@ kTF1_red
Uses TF1_red for interpolation.
@ kTSpline3_red
Uses TSpline3_red for interpolation.
TResponseFunction_red * CreateResponseFunction(TGraph *&graph, const RespFuncType SplineRespFuncType, const SplineInterpolation SplineInterpolationType, const std::string &Title)
KS: Create Response Function using TGraph.
void ApplyKnotWeightCap(TGraph *xsecgraph, const int splineParsIndex, ParameterHandlerGeneric *XsecCov)
EM: Apply capping to knot weight for specified spline parameter. param graph needs to have been set i...
void ApplyKnotWeightCapTSpline3(TSpline3 *&Spline, const int splineParsIndex, ParameterHandlerGeneric *XsecCov)
EM: Apply capping to knot weight for specified spline parameter. param graph needs to have been set i...
std::vector< std::vector< TF1_red * > > ReduceTF1(std::vector< std::vector< TF1 * > > &MasterSpline)
CW: Reduced the TF1 to TF1_red.
std::vector< std::vector< TSpline3_red * > > ReduceTSpline3(std::vector< std::vector< TSpline3 * > > &MasterSpline)
CW: Reduced the TSpline3 to TSpline3_red.
bool isFlat(TSpline3_red *&spl)
CW: Helper function used in the constructor, tests to see if the spline is flat.
Custom exception class used throughout MaCh3.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
Class responsible for handling of systematic error parameters with different types defined in the con...
double GetParSplineKnotUpperBound(const int i) const
EM: value at which we cap spline knot weight.
double GetParSplineKnotLowerBound(const int i) const
EM: value at which we cap spline knot weight.
CW: A reduced TF1 class only. Only saves parameters for each TF1 and how many parameters each paramet...
void SetSize(M3::int_t nSpline)
Set the size.
void SetFunc(TF1 *&Func)
Set the function.
TF1_red()
Empty constructor.
void Print() override
Print detailed info.
M3::int_t GetNp() override
DL: Get number of points.
TF1 * ConstructTF1(const std::string &function, const int xmin, const int xmax)
KS: Make a TF1 from the reduced TF1.
TF1_red(TF1 *&Function)
The TF1 constructor with deep copy.
M3::float_t * Par
The parameters.
virtual ~TF1_red()
Empty destructor.
double Eval(const double var) override
Evaluate a variation.
TF1_red(M3::int_t nSize, M3::float_t *Array)
The useful constructor with deep copy.
double GetParameter(M3::int_t Parameter)
Get a parameter value.
void SetParameter(M3::int_t Parameter, M3::float_t Value)
Set a parameter to a value.
int GetSize()
Get the size.
KS: A reduced ResponseFunction Generic function used for evaluating weight.
virtual ~TResponseFunction_red()
Empty destructor.
TResponseFunction_red()
Empty constructor.
virtual M3::int_t GetNp()=0
DL: Get number of points.
virtual double Eval(const double var)=0
Evaluate a variation.
virtual void Print()=0
KS: Printer.
CW: Reduced TSpline3 class.
double Eval(double var) override
CW: Evaluate the weight from a variation.
void GetKnot(int i, M3::float_t &xtmp, M3::float_t &ytmp)
int FindX(double x)
Find the segment relevant to this variation in x.
virtual ~TSpline3_red()
Empty destructor.
M3::float_t ** Par
Always uses a third order polynomial, so hard-code the number of coefficients in implementation.
TSpline3_red(M3::float_t *X, M3::float_t *Y, M3::int_t N, M3::float_t **P)
constructor taking parameters
M3::int_t nPoints
Number of points/knot in TSpline3.
TSpline3_red(TSpline3 *&spline, SplineInterpolation InterPolation=kTSpline3)
The constructor that takes a TSpline3 pointer and copies in to memory.
M3::int_t GetNp() override
CW: Get the number of points.
TSpline3 * ConstructTSpline3()
CW: Make a TSpline3 from the reduced splines.
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.
void Print() override
Print detailed info.
void SetFunc(TSpline3 *&spline, SplineInterpolation InterPolation=kTSpline3)
Set the function .
M3::float_t * XPos
Positions of each x for each knot.
TSpline3_red()
Empty constructor.
M3::float_t * YResp
y-value for each knot
constexpr static const double DefSplineKnotUpBound
Default value for spline knot capping, default mean not capping is being applied.
constexpr static const double DefSplineKnotLowBound
Default value for spline knot capping, default mean not capping is being applied.
CW: Add a struct to hold info about the splinified xsec parameters and help with FindSplineSegment.
virtual ~FastSplineInfo()=default
Destructor.
std::vector< M3::float_t > xPts
Array of the knots positions.
const double * splineParsPointer
Array of the knots positions.
M3::int_t CurrSegment
Array of what segment of spline we're currently interested in. Gets updated once per MCMC iteration.
FastSplineInfo()
Constructor.
M3::int_t nPts
Number of points in spline.