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"
38 std::vector<M3::float_t>
xPts;
52 if(xsecgraph ==
nullptr){
53 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));
62 for(
int knotId = 0; knotId < xsecgraph->GetN(); knotId++){
66 if( xsecgraph->GetPoint(knotId, x, y) == -1) {
74 xsecgraph->SetPoint(knotId, x, y);
79 for(
int knotId = 0; knotId < xsecgraph->GetN(); knotId++){
83 if( xsecgraph->GetPoint(knotId, x, y) == -1) {
87 if(std::abs(y - 1.0) > 1e-5)
isFlat =
false;
92 xsecgraph->SetPoint(0, 0.0, 1.0);
101 if(Spline ==
nullptr) {
102 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));
106 std::string oldName = Spline->GetName();
110 const int NValues = Spline->GetNp();
111 std::vector<double> XVals(NValues);
112 std::vector<double> YVals(NValues);
113 for(
int knotId = 0; knotId < NValues; knotId++){
116 Spline->GetKnot(knotId, x, y);
126 Spline =
new TSpline3(oldName.c_str(), XVals.data(), YVals.data(), NValues);
140 virtual double Eval(
const double var)=0;
169 for (
int i = 0; i <
length; ++i) {
183 if (
Par !=
nullptr)
delete[]
Par;
185 for (
int i = 0; i <
length; ++i) {
193 inline double Eval(
const double var)
override {
217 Par[Parameter] = Value;
223 MACH3LOG_ERROR(
"You requested parameter number {} but length is {} parameters", Parameter,
length);
227 return Par[Parameter];
241 for (
int i = 0; i <
length; i++) {
247 inline TF1*
ConstructTF1(
const std::string&
function,
const int xmin,
const int xmax) {
248 TF1 *func =
new TF1(
"TF1",
function.c_str(), xmin, xmax);
249 for(
int i = 0; i <
length; ++i) {
250 func->SetParameter(i,
Par[i]);
282 SetFunc(spline, InterPolation);
294 for(
int j = 0; j < N; ++j){
302 if((
Par[j][0] == -999) | (
Par[j][1] ==-999) | (
Par[j][2] ==-999) | (
XPos[j] ==-999) | (
YResp[j] ==-999)){
303 MACH3LOG_ERROR(
"******************* Bad parameter values when constructing TSpline3_red *********************");
304 MACH3LOG_ERROR(
"Passed values (i, x, y, b, c, d): {}, {}, {}, {}, {}, {}", j, X[j], Y[j], P[j][0], P[j][1], P[j][2]);
305 MACH3LOG_ERROR(
"Set values (i, x, y, b, c, d): {}, {}, {}, {}, {}, {}", j,
XPos[j],
YResp[j],
Par[j][0],
Par[j][1],
Par[j][2]);
306 MACH3LOG_ERROR(
"*********************************************************************************************");
314 if (
Par !=
nullptr) {
315 for (
int i = 0; i <
nPoints; ++i) {
334 for (
int i = 0; i <
nPoints; ++i) {
337 double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
338 spline->GetCoeff(i, x, y, b, c, d);
354 for (
int k = 0; k <
nPoints; ++k) {
357 Double_t x1, y1, b1, c1, d1, x2, y2, b2, c2, d2 = 0;
358 spline->GetCoeff(k, x1, y1, b1, c1, d1);
364 spline->GetCoeff(k + 1, x2, y2, b2, c2, d2);
365 tempb = (y2-y1)/(x2-x1);
376 else if(InterPolation ==
kAkima)
379 for (
int i = 0; i <
nPoints; ++i) {
383 double x = -999.99, y = -999.99;
384 spline->GetKnot(i, x, y);
393 for (
int i = -2; i <=
nPoints; ++i) {
414 for(
int i =2; i<=
nPoints+2; i++){
415 if (std::abs(mvals[i+1] - mvals[i]) + std::abs(mvals[i-1] - mvals[i-2]) != 0.0){
416 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]));
418 else{svals[i-2] = mvals[i];}
422 for(
int i = 0; i <
nPoints; i++){
436 for(
int i = 0; i <
nPoints-1; i++){
437 double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
438 spline->GetCoeff(i, x, y, b, c, d);
440 if((c == 0.0 && d == 0.0)){
459 for (
int i = 0; i <
nPoints; ++i) {
463 double x = -999.99, y = -999.99;
464 spline->GetKnot(i, x, y);
486 for (
int i = 0; i <
nPoints-1; ++i) {
491 Tangents[0] = Secants[0];
498 for (
int i = 1; i <
nPoints-1; ++i) {
499 if ((Secants[i-1] >= 0.0 && Secants[i] >= 0.0) | (Secants[i-1] < 0.0 && Secants[i] < 0.0)){
500 Tangents[i] =
M3::float_t((Secants[i-1] + Secants[i]) /2.0);
505 for (
int i = 0; i <
nPoints-1; ++i) {
506 if (Secants[i] == 0.0){
512 alpha = Tangents[i] / Secants[i];
513 beta = Tangents[i+1] / Secants[i];
522 if (alpha * alpha + beta * beta >9.0){
524 Tangents[i] = tau * alpha * Secants[i];
525 Tangents[i+1] = tau * beta * Secants[i];
531 for(
int i = 0; i <
nPoints-1; i++){
535 b = Tangents[i] * dx;
540 Par[i][1] = c / (dx * dx);
541 Par[i][2] = d / (dx * dx * dx);
543 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)){
544 MACH3LOG_INFO(
"Bad spline parameters for segment {}: (b, c, d) = {}, {}, {}. This will cause problems with GPU.",
560 for(
int i = 0; i <
nPoints-1; i++){
561 double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
562 spline->GetCoeff(i, x, y, b, c, d);
564 if((c == 0.0 && d == 0.0)){
578 for (
int i = 0; i <
nPoints; ++i)
581 double x = -999.99, y = -999.99;
582 spline->GetKnot(i, x, y);
608 for (
int i = 1; i <
nPoints - 1; ++i)
613 M3::float_t term1 = (1.0 - T) * (1.0 + B) * (1.0 + C) * 0.5 * d0;
614 M3::float_t term2 = (1.0 - T) * (1.0 - B) * (1.0 - C) * 0.5 * d1;
616 Tangents[i] = term1 + term2;
624 for (
int i = 0; i <
nPoints - 1; ++i)
632 Par[i][1] = (3*dy/(dx*dx)) - (2*m0 + m1)/dx;
633 Par[i][2] = (m0 + m1 - 2*dy/dx) / (dx*dx);
636 i, dx, dy,
Par[i][0],
Par[i][1],
Par[i][2]);
647 MACH3LOG_ERROR(
"Unsupported interpolation type {}",
static_cast<int>(InterPolation));
658 for (
int i = 0; i <
nPoints; ++i) {
659 if (
Par[i] !=
nullptr) {
691 while (kHigh - segment > 1) {
693 kHalf = (segment + kHigh)/2;
695 if (x >
XPos[kHalf]) {
708 inline double Eval(
double var)
override {
710 int segment =
FindX(var);
716 double weight = y+dx*(b+dx*(c+d*dx));
740 #ifdef _LOW_MEMORY_STRUCTS_
741 std::vector<Double_t> xPosDoubles(
nPoints);
742 std::vector<Double_t> yPosDoubles(
nPoints);
743 for (Int_t i = 0; i <
nPoints; ++i) {
744 xPosDoubles[i] =
static_cast<Double_t
>(
XPos[i]);
745 yPosDoubles[i] =
static_cast<Double_t
>(
YResp[i]);
747 TSpline3 *spline =
new TSpline3(
"Spline", xPosDoubles.data(), yPosDoubles.data(),
static_cast<int>(
nPoints));
751 for (Int_t i = 0; i <
nPoints; ++i) {
752 spline->SetPointCoeff(i,
Par[i][0],
Par[i][1],
Par[i][2]);
762 for (
int i = 0; i <
nPoints; ++i) {
784 int Np = spl->
GetNp();
789 for(
int i = 0; i < Np; i++) {
801 inline std::vector<std::vector<TSpline3_red*> >
ReduceTSpline3(std::vector<std::vector<TSpline3*> > &MasterSpline) {
803 std::vector<std::vector<TSpline3*> >::iterator OuterIt;
804 std::vector<TSpline3*>::iterator InnerIt;
807 std::vector<std::vector<TSpline3_red*> > ReducedVector;
808 ReducedVector.reserve(MasterSpline.size());
811 int OuterCounter = 0;
812 for (OuterIt = MasterSpline.begin(); OuterIt != MasterSpline.end(); ++OuterIt, ++OuterCounter) {
814 std::vector<TSpline3_red*> TempVector;
815 TempVector.reserve(OuterIt->size());
816 int InnerCounter = 0;
818 for (InnerIt = OuterIt->begin(); InnerIt != OuterIt->end(); ++InnerIt, ++InnerCounter) {
820 TSpline3 *spline = (*InnerIt);
823 if (spline != NULL) {
828 TempVector.push_back(red);
830 ReducedVector.push_back(TempVector);
833 return ReducedVector;
839 inline std::vector<std::vector<TF1_red*> >
ReduceTF1(std::vector<std::vector<TF1*> > &MasterSpline) {
841 std::vector<std::vector<TF1*> >::iterator OuterIt;
842 std::vector<TF1*>::iterator InnerIt;
845 std::vector<std::vector<TF1_red*> > ReducedVector;
846 ReducedVector.reserve(MasterSpline.size());
849 int OuterCounter = 0;
850 for (OuterIt = MasterSpline.begin(); OuterIt != MasterSpline.end(); ++OuterIt, ++OuterCounter) {
852 std::vector<TF1_red*> TempVector;
853 TempVector.reserve(OuterIt->size());
854 int InnerCounter = 0;
856 for (InnerIt = OuterIt->begin(); InnerIt != OuterIt->end(); ++InnerIt, ++InnerCounter) {
858 TF1* spline = (*InnerIt);
861 if (spline != NULL) {
866 TempVector.push_back(red);
868 ReducedVector.push_back(TempVector);
871 return ReducedVector;
883 const std::string& Title) {
887 if (graph && graph->GetN() > 1)
892 TSpline3* spline =
nullptr;
896 spline =
new TSpline3(Title.c_str(), graph);
897 spline->SetNameTitle(Title.c_str(), Title.c_str());
900 spline_red =
new TSpline3_red(spline, SplineInterpolationType);
902 RespFunc = spline_red;
904 else if(SplineRespFuncType ==
kTF1_red)
907 TF1 *Fitter =
nullptr;
910 if(graph->GetN() != 2) {
911 MACH3LOG_ERROR(
"Trying to make TF1 from more than 2 knots. Knots = {}", graph->GetN());
917 Fitter =
new TF1(Title.c_str(),
"([1]+[0]*x)", graph->GetX()[0], graph->GetX()[graph->GetN()-1]);
926 graph->Fit(Fitter,
"Q0");
945 #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.
const M3::float_t * splineParsPointer
Array of the knots positions.
std::vector< M3::float_t > xPts
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.